Fork me on GitHub

source: git/display/DelphesPlotSummary.cc@ 53b78e8

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 53b78e8 was 53b78e8, checked in by Christophe Delaere <christophe.delaere@…>, 10 years ago

New signal to show the progress of plot processing

Progress is now shown in the status bar when generating std plots.

  • Property mode set to 100644
File size: 6.4 KB
RevLine 
[2ca23b5]1#include "display/DelphesPlotSummary.h"
[3f51314]2#include "TRootEmbeddedCanvas.h"
[ad015db]3#include <algorithm>
[2ca23b5]4
[ad015db]5bool vecsorter (TLorentzVector i,TLorentzVector j) { return (i.Pt()>j.Pt()); }
[3f51314]6
7DelphesPlotSummary::DelphesPlotSummary(TEveWindowTab* tab):tab_(tab) {}
[2ca23b5]8
9DelphesPlotSummary::~DelphesPlotSummary() {}
10
[53b78e8]11void DelphesPlotSummary::Progress(Int_t p)
12{
13 Emit("Progress(Int_t)",p);
14}
15
[3f51314]16void DelphesPlotSummary::Init(std::vector<DelphesBranchBase*>& elements) {
17 elements_ = &elements;
18 // loop on the elements, and create tabs
19 for(std::vector<DelphesBranchBase*>::iterator data=elements.begin();data<elements.end();++data) {
20 // the canvas
21 TEveWindowSlot* slot = tab_->NewSlot();
22 TRootEmbeddedCanvas* trec = new TRootEmbeddedCanvas();
23 TCanvas* canvas = trec->GetCanvas();
24 TEveWindowFrame * wf = slot->MakeFrame(trec);
25 wf->SetElementName((*data)->GetName());
26 canvas->Divide(3,3);
27 canvases_[(*data)->GetName()] = canvas;
28 // the histograms
29 TH1F* h;
30 std::vector<TH1F*> histograms;
31 histograms.reserve(9);
32 h = new TH1F(Form("%sPt",(*data)->GetName()),Form("%s Pt",(*data)->GetName()),100,0,-1);
33 histograms.push_back(h);
34 h = new TH1F(Form("%sEta",(*data)->GetName()),Form("%s Eta",(*data)->GetName()),100,0,-1);
35 histograms.push_back(h);
36 h = new TH1F(Form("%sPhi",(*data)->GetName()),Form("%s Phi",(*data)->GetName()),100,0,-1);
37 histograms.push_back(h);
38 h = new TH1F(Form("l%sPt",(*data)->GetName()),Form("leading %s Pt",(*data)->GetName()),100,0,-1);
39 histograms.push_back(h);
40 h = new TH1F(Form("l%sEta",(*data)->GetName()),Form("leading %s Eta",(*data)->GetName()),100,0,-1);
41 histograms.push_back(h);
42 h = new TH1F(Form("l%sPhi",(*data)->GetName()),Form("leading %s Phi",(*data)->GetName()),100,0,-1);
43 histograms.push_back(h);
44 h = new TH1F(Form("sl%sPt",(*data)->GetName()),Form("subleading %s Pt",(*data)->GetName()),100,0,-1);
45 histograms.push_back(h);
46 h = new TH1F(Form("sl%sEta",(*data)->GetName()),Form("subleading %s Eta",(*data)->GetName()),100,0,-1);
47 histograms.push_back(h);
48 h = new TH1F(Form("sl%sPhi",(*data)->GetName()),Form("subleading %s Phi",(*data)->GetName()),100,0,-1);
49 histograms.push_back(h);
50 histograms_[(*data)->GetName()] = histograms;
51 }
52}
53
54void DelphesPlotSummary::FillSample(ExRootTreeReader* treeReader, Int_t event_id) {
[53b78e8]55 Int_t entries = treeReader->GetEntries();
56 for(Int_t i=0;i<entries;++i) {
[3f51314]57 treeReader->ReadEntry(i);
58 for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
59 std::vector<TLorentzVector> vectors = (*element)->GetVectors();
[ad015db]60 std::sort(vectors.begin(), vectors.end(), vecsorter);
[3f51314]61 std::vector<TH1F*> histograms = histograms_[(*element)->GetName()];
62 for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
63 histograms[0]->Fill(it->Pt());
64 histograms[1]->Fill(it->Eta());
65 histograms[2]->Fill(it->Phi());
66 if(it==vectors.begin()) {
67 histograms[3]->Fill(it->Pt());
68 histograms[4]->Fill(it->Eta());
69 histograms[5]->Fill(it->Phi());
70 }
71 if(it==vectors.begin()+1) {
72 histograms[6]->Fill(it->Pt());
73 histograms[7]->Fill(it->Eta());
74 histograms[8]->Fill(it->Phi());
75 }
76 }
77 }
[53b78e8]78 Progress(int(100*i/entries));
[3f51314]79 }
80 treeReader->ReadEntry(event_id);
[53b78e8]81 Progress(100);
[3f51314]82}
[2ca23b5]83
[3f51314]84void DelphesPlotSummary::Draw() {
85 for(std::map< TString, TCanvas* >::iterator it=canvases_.begin(); it!=canvases_.end(); ++it) {
86 TCanvas* c = it->second;
87 std::vector<TH1F*> histograms = histograms_[it->first];
88 std::vector<TH1F*> eventProfiles = eventProfiles_[it->first];
89 std::vector<TMarker*> eventMarkers = eventMarkers_[it->first];
90 for(Int_t i=0;i<9;++i) {
91 c->cd(i+1);
[2695ae1]92 if(histograms[i]->GetEntries()==0) continue;
[3f51314]93 histograms[i]->Draw();
[ad015db]94 if(i<3) {
[3f51314]95 eventProfiles[i]->Draw("same");
[ad015db]96 } else {
97 if(eventMarkers.size()>(unsigned int)(i-3)) {
98 eventMarkers[i-3]->Draw("same");
99 }
100 }
[3f51314]101 }
[ad015db]102 c->Update();
[3f51314]103 }
104}
[2ca23b5]105
[6301e02]106void DelphesPlotSummary::FillEvent() { //TODO make it faster... try to reuse objects and simply change the values
[3f51314]107 // clear previous markers
108 for(std::map< TString, std::vector<TMarker*> >::iterator mv = eventMarkers_.begin(); mv!=eventMarkers_.end();++mv) {
109 for(std::vector<TMarker*>::iterator m = mv->second.begin(); m<mv->second.end();++m) {
110 delete *m;
111 }
112 }
113 eventMarkers_.clear();
114 for(std::map< TString, std::vector<TH1F*> >::iterator hv = eventProfiles_.begin(); hv!=eventProfiles_.end();++hv) {
115 for(std::vector<TH1F*>::iterator h = hv->second.begin(); h<hv->second.end();++h) {
116 delete *h;
117 }
118 }
119 eventProfiles_.clear();
120 // loop over the elements and fill markers with event data
121 TMarker *m;
122 for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
[ad015db]123 std::vector<TMarker*> mv;
[3f51314]124 std::vector<TLorentzVector> vectors = (*element)->GetVectors();
[ad015db]125 std::sort(vectors.begin(), vectors.end(), vecsorter);
[3f51314]126 std::vector<TH1F*> histograms = histograms_[(*element)->GetName()];
127 TH1F* h1 = (TH1F*)histograms[0]->Clone(); h1->Reset(); h1->SetLineColor(kBlue);
128 TH1F* h2 = (TH1F*)histograms[1]->Clone(); h2->Reset(); h2->SetLineColor(kBlue);
129 TH1F* h3 = (TH1F*)histograms[2]->Clone(); h3->Reset(); h3->SetLineColor(kBlue);
130 for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
131 h1->Fill(it->Pt());
132 h2->Fill(it->Eta());
133 h3->Fill(it->Phi());
134 if(it==vectors.begin()) {
[ad015db]135 m = new TMarker(it->Pt(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]136 mv.push_back(m);
[ad015db]137 m = new TMarker(it->Eta(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]138 mv.push_back(m);
[ad015db]139 m = new TMarker(it->Phi(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]140 mv.push_back(m);
141 }
142 if(it==vectors.begin()+1) {
[ad015db]143 m = new TMarker(it->Pt(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]144 mv.push_back(m);
[ad015db]145 m = new TMarker(it->Eta(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]146 mv.push_back(m);
[ad015db]147 m = new TMarker(it->Phi(),0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
[3f51314]148 mv.push_back(m);
149 }
150 }
151 std::vector<TH1F*> hv;
152 hv.push_back(h1);
153 hv.push_back(h2);
154 hv.push_back(h3);
155 eventProfiles_[(*element)->GetName()] = hv;
156 eventMarkers_[(*element)->GetName()] = mv;
157 }
158}
Note: See TracBrowser for help on using the repository browser.