Fork me on GitHub

source: git/display/DelphesPlotSummary.cc@ 6301e02

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

Better GUI

Improved the event navigation widget + batch operations

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