Fork me on GitHub

source: git/display/DelphesPlotSummary.cc@ 2695ae1

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

Avoid error messages for empty plots

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