Fork me on GitHub

Changeset 3f51314 in git for display/DelphesPlotSummary.cc


Ignore:
Timestamp:
Oct 25, 2014, 11:32:13 PM (10 years ago)
Author:
Christophe Delaere <christophe.delaere@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
ad015db
Parents:
2ca23b5
Message:

First incomplete version with summary plots

Most of the functionalities are there, but it is not yet properly
integrated in the GUI.
Known issues:

  • size of the event markers
  • all markers appear the same (missing index?)
File:
1 edited

Legend:

Unmodified
Added
Removed
  • display/DelphesPlotSummary.cc

    r2ca23b5 r3f51314  
    11#include "display/DelphesPlotSummary.h"
     2#include "TRootEmbeddedCanvas.h"
    23
    3 DelphesPlotSummary::DelphesPlotSummary(TEveWindowTab* tab) {}
     4
     5DelphesPlotSummary::DelphesPlotSummary(TEveWindowTab* tab):tab_(tab) {}
    46
    57DelphesPlotSummary::~DelphesPlotSummary() {}
    68
    7 void DelphesPlotSummary::Init(std::vector<DelphesBranchBase*>& elements) {}
     9void DelphesPlotSummary::Init(std::vector<DelphesBranchBase*>& elements) {
     10  elements_ = &elements;
     11  // loop on the elements, and create tabs
     12  for(std::vector<DelphesBranchBase*>::iterator data=elements.begin();data<elements.end();++data) {
     13    // the canvas
     14    TEveWindowSlot* slot = tab_->NewSlot();
     15    TRootEmbeddedCanvas* trec = new TRootEmbeddedCanvas();
     16    TCanvas* canvas = trec->GetCanvas();
     17    TEveWindowFrame * wf = slot->MakeFrame(trec);
     18    wf->SetElementName((*data)->GetName());
     19    canvas->Divide(3,3);
     20    canvases_[(*data)->GetName()] = canvas;
     21    // the histograms
     22    TH1F* h;
     23    std::vector<TH1F*> histograms;
     24    histograms.reserve(9);
     25    h = new TH1F(Form("%sPt",(*data)->GetName()),Form("%s Pt",(*data)->GetName()),100,0,-1);
     26    histograms.push_back(h);
     27    h = new TH1F(Form("%sEta",(*data)->GetName()),Form("%s Eta",(*data)->GetName()),100,0,-1);
     28    histograms.push_back(h);
     29    h = new TH1F(Form("%sPhi",(*data)->GetName()),Form("%s Phi",(*data)->GetName()),100,0,-1);
     30    histograms.push_back(h);
     31    h = new TH1F(Form("l%sPt",(*data)->GetName()),Form("leading %s Pt",(*data)->GetName()),100,0,-1);
     32    histograms.push_back(h);
     33    h = new TH1F(Form("l%sEta",(*data)->GetName()),Form("leading %s Eta",(*data)->GetName()),100,0,-1);
     34    histograms.push_back(h);
     35    h = new TH1F(Form("l%sPhi",(*data)->GetName()),Form("leading %s Phi",(*data)->GetName()),100,0,-1);
     36    histograms.push_back(h);
     37    h = new TH1F(Form("sl%sPt",(*data)->GetName()),Form("subleading %s Pt",(*data)->GetName()),100,0,-1);
     38    histograms.push_back(h);
     39    h = new TH1F(Form("sl%sEta",(*data)->GetName()),Form("subleading %s Eta",(*data)->GetName()),100,0,-1);
     40    histograms.push_back(h);
     41    h = new TH1F(Form("sl%sPhi",(*data)->GetName()),Form("subleading %s Phi",(*data)->GetName()),100,0,-1);
     42    histograms.push_back(h);
     43    histograms_[(*data)->GetName()] = histograms;
     44  }
     45}
    846
    9 void DelphesPlotSummary::FillSample(ExRootTreeReader* treeReader) {}
     47void DelphesPlotSummary::FillSample(ExRootTreeReader* treeReader, Int_t event_id) {
     48  for(Int_t i=0;i<treeReader->GetEntries();++i) {
     49    treeReader->ReadEntry(i);
     50    for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
     51      std::vector<TLorentzVector> vectors = (*element)->GetVectors();
     52      std::vector<TH1F*> histograms = histograms_[(*element)->GetName()];
     53      for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
     54        histograms[0]->Fill(it->Pt());
     55        histograms[1]->Fill(it->Eta());
     56        histograms[2]->Fill(it->Phi());
     57        if(it==vectors.begin()) {
     58          histograms[3]->Fill(it->Pt());
     59          histograms[4]->Fill(it->Eta());
     60          histograms[5]->Fill(it->Phi());
     61        }
     62        if(it==vectors.begin()+1) {
     63          histograms[6]->Fill(it->Pt());
     64          histograms[7]->Fill(it->Eta());
     65          histograms[8]->Fill(it->Phi());
     66        }
     67      }
     68    }
     69  }
     70  treeReader->ReadEntry(event_id);
     71}
    1072
    11 void DelphesPlotSummary::FillEvent() {}
     73void DelphesPlotSummary::Draw() {
     74  for(std::map< TString, TCanvas* >::iterator it=canvases_.begin(); it!=canvases_.end(); ++it) {
     75    TCanvas* c = it->second;
     76    std::vector<TH1F*> histograms = histograms_[it->first];
     77    std::vector<TH1F*> eventProfiles = eventProfiles_[it->first];
     78    std::vector<TMarker*> eventMarkers = eventMarkers_[it->first];
     79    for(Int_t i=0;i<9;++i) {
     80      c->cd(i+1);
     81      histograms[i]->Draw();
     82      if(i<3)
     83        eventProfiles[i]->Draw("same");
     84      else
     85        eventMarkers[i-3]->Draw("same");
     86    }
     87  }
     88}
     89
     90void DelphesPlotSummary::FillEvent() {
     91  // clear previous markers
     92  for(std::map< TString, std::vector<TMarker*> >::iterator mv = eventMarkers_.begin(); mv!=eventMarkers_.end();++mv) {
     93    for(std::vector<TMarker*>::iterator m = mv->second.begin(); m<mv->second.end();++m) {
     94      delete *m;
     95    }
     96  }
     97  eventMarkers_.clear();
     98  for(std::map< TString, std::vector<TH1F*> >::iterator hv = eventProfiles_.begin(); hv!=eventProfiles_.end();++hv) {
     99    for(std::vector<TH1F*>::iterator h = hv->second.begin(); h<hv->second.end();++h) {
     100      delete *h;
     101    }
     102  }
     103  eventProfiles_.clear();
     104  // loop over the elements and fill markers with event data
     105  TMarker *m;
     106  std::vector<TMarker*> mv;
     107  for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
     108    std::vector<TLorentzVector> vectors = (*element)->GetVectors();
     109    std::vector<TH1F*> histograms = histograms_[(*element)->GetName()];
     110    TH1F* h1 = (TH1F*)histograms[0]->Clone(); h1->Reset(); h1->SetLineColor(kBlue);
     111    TH1F* h2 = (TH1F*)histograms[1]->Clone(); h2->Reset(); h2->SetLineColor(kBlue);
     112    TH1F* h3 = (TH1F*)histograms[2]->Clone(); h3->Reset(); h3->SetLineColor(kBlue);
     113    for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
     114      h1->Fill(it->Pt());
     115      h2->Fill(it->Eta());
     116      h3->Fill(it->Phi());
     117      if(it==vectors.begin()) {
     118        m = new TMarker(it->Pt(),0,29); m->SetMarkerColor(kBlue);
     119        mv.push_back(m);
     120        m = new TMarker(it->Eta(),0,29); m->SetMarkerColor(kBlue);
     121        mv.push_back(m);
     122        m = new TMarker(it->Phi(),0,29); m->SetMarkerColor(kBlue);
     123        mv.push_back(m);
     124      }
     125      if(it==vectors.begin()+1) {
     126        m = new TMarker(it->Pt(),0,29); m->SetMarkerColor(kBlue);
     127        mv.push_back(m);
     128        m = new TMarker(it->Eta(),0,29); m->SetMarkerColor(kBlue);
     129        mv.push_back(m);
     130        m = new TMarker(it->Phi(),0,29); m->SetMarkerColor(kBlue);
     131        mv.push_back(m);
     132      }
     133    }
     134    std::vector<TH1F*> hv;
     135    hv.push_back(h1);
     136    hv.push_back(h2);
     137    hv.push_back(h3);
     138    eventProfiles_[(*element)->GetName()] = hv;
     139    eventMarkers_[(*element)->GetName()] = mv;
     140  }
     141}
Note: See TracChangeset for help on using the changeset viewer.