1 | /*
|
---|
2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
4 | *
|
---|
5 | * This program is free software: you can redistribute it and/or modify
|
---|
6 | * it under the terms of the GNU General Public License as published by
|
---|
7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
8 | * (at your option) any later version.
|
---|
9 | *
|
---|
10 | * This program is distributed in the hope that it will be useful,
|
---|
11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | * GNU General Public License for more details.
|
---|
14 | *
|
---|
15 | * You should have received a copy of the GNU General Public License
|
---|
16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
17 | */
|
---|
18 |
|
---|
19 | #include "display/DelphesPlotSummary.h"
|
---|
20 | #include "TRootEmbeddedCanvas.h"
|
---|
21 | #include <algorithm>
|
---|
22 |
|
---|
23 | bool vecsorter (TLorentzVector i,TLorentzVector j) { return (i.Pt()>j.Pt()); }
|
---|
24 |
|
---|
25 | DelphesPlotSummary::DelphesPlotSummary(TEveWindowTab* tab):tab_(tab) {}
|
---|
26 |
|
---|
27 | DelphesPlotSummary::~DelphesPlotSummary() {}
|
---|
28 |
|
---|
29 | void DelphesPlotSummary::Progress(Int_t p)
|
---|
30 | {
|
---|
31 | Emit("Progress(Int_t)",p);
|
---|
32 | }
|
---|
33 |
|
---|
34 | void DelphesPlotSummary::Init(std::vector<DelphesBranchBase*>& elements) {
|
---|
35 | elements_ = &elements;
|
---|
36 | // loop on the elements, and create tabs
|
---|
37 | for(std::vector<DelphesBranchBase*>::iterator data=elements.begin();data<elements.end();++data) {
|
---|
38 | // the canvas
|
---|
39 | TEveWindowSlot* slot = tab_->NewSlot();
|
---|
40 | TRootEmbeddedCanvas* trec = new TRootEmbeddedCanvas();
|
---|
41 | TCanvas* canvas = trec->GetCanvas();
|
---|
42 | TEveWindowFrame * wf = slot->MakeFrame(trec);
|
---|
43 | wf->SetElementName((*data)->GetName());
|
---|
44 | canvas->Divide(3,3);
|
---|
45 | canvases_[(*data)->GetName()] = canvas;
|
---|
46 | // the histograms
|
---|
47 | TH1F* h;
|
---|
48 | std::vector<TH1F*> histograms;
|
---|
49 | histograms.reserve(9);
|
---|
50 | h = new TH1F(Form("%sPt",(*data)->GetName()),Form("%s Pt",(*data)->GetName()),100,0,-1);
|
---|
51 | histograms.push_back(h);
|
---|
52 | h = new TH1F(Form("%sEta",(*data)->GetName()),Form("%s Eta",(*data)->GetName()),100,0,-1);
|
---|
53 | histograms.push_back(h);
|
---|
54 | h = new TH1F(Form("%sPhi",(*data)->GetName()),Form("%s Phi",(*data)->GetName()),100,0,-1);
|
---|
55 | histograms.push_back(h);
|
---|
56 | h = new TH1F(Form("l%sPt",(*data)->GetName()),Form("leading %s Pt",(*data)->GetName()),100,0,-1);
|
---|
57 | histograms.push_back(h);
|
---|
58 | h = new TH1F(Form("l%sEta",(*data)->GetName()),Form("leading %s Eta",(*data)->GetName()),100,0,-1);
|
---|
59 | histograms.push_back(h);
|
---|
60 | h = new TH1F(Form("l%sPhi",(*data)->GetName()),Form("leading %s Phi",(*data)->GetName()),100,0,-1);
|
---|
61 | histograms.push_back(h);
|
---|
62 | h = new TH1F(Form("sl%sPt",(*data)->GetName()),Form("subleading %s Pt",(*data)->GetName()),100,0,-1);
|
---|
63 | histograms.push_back(h);
|
---|
64 | h = new TH1F(Form("sl%sEta",(*data)->GetName()),Form("subleading %s Eta",(*data)->GetName()),100,0,-1);
|
---|
65 | histograms.push_back(h);
|
---|
66 | h = new TH1F(Form("sl%sPhi",(*data)->GetName()),Form("subleading %s Phi",(*data)->GetName()),100,0,-1);
|
---|
67 | histograms.push_back(h);
|
---|
68 | histograms_[(*data)->GetName()] = histograms;
|
---|
69 | // the event histograms
|
---|
70 | TH1F* h1 = (TH1F*)histograms[0]->Clone(); h1->Reset(); h1->SetLineColor(kBlue);
|
---|
71 | TH1F* h2 = (TH1F*)histograms[1]->Clone(); h2->Reset(); h2->SetLineColor(kBlue);
|
---|
72 | TH1F* h3 = (TH1F*)histograms[2]->Clone(); h3->Reset(); h3->SetLineColor(kBlue);
|
---|
73 | std::vector<TH1F*> hv;
|
---|
74 | hv.push_back(h1);
|
---|
75 | hv.push_back(h2);
|
---|
76 | hv.push_back(h3);
|
---|
77 | eventProfiles_[(*data)->GetName()] = hv;
|
---|
78 | // the event markers
|
---|
79 | TMarker *m;
|
---|
80 | std::vector<TMarker*> mv;
|
---|
81 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
82 | mv.push_back(m);
|
---|
83 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
84 | mv.push_back(m);
|
---|
85 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
86 | mv.push_back(m);
|
---|
87 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
88 | mv.push_back(m);
|
---|
89 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
90 | mv.push_back(m);
|
---|
91 | m = new TMarker(0,0,29); m->SetMarkerColor(kBlue); m->SetMarkerSize(3);
|
---|
92 | mv.push_back(m);
|
---|
93 | eventMarkers_[(*data)->GetName()] = mv;
|
---|
94 | }
|
---|
95 | }
|
---|
96 |
|
---|
97 | void DelphesPlotSummary::FillSample(ExRootTreeReader* treeReader, Int_t event_id) {
|
---|
98 | Int_t entries = treeReader->GetEntries();
|
---|
99 | for(Int_t i=0;i<entries;++i) {
|
---|
100 | treeReader->ReadEntry(i);
|
---|
101 | for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
|
---|
102 | std::vector<TLorentzVector> vectors = (*element)->GetVectors();
|
---|
103 | std::sort(vectors.begin(), vectors.end(), vecsorter);
|
---|
104 | std::vector<TH1F*> histograms = histograms_[(*element)->GetName()];
|
---|
105 | for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
|
---|
106 | histograms[0]->Fill(it->Pt());
|
---|
107 | histograms[1]->Fill(it->Eta());
|
---|
108 | histograms[2]->Fill(it->Phi());
|
---|
109 | if(it==vectors.begin()) {
|
---|
110 | histograms[3]->Fill(it->Pt());
|
---|
111 | histograms[4]->Fill(it->Eta());
|
---|
112 | histograms[5]->Fill(it->Phi());
|
---|
113 | }
|
---|
114 | if(it==vectors.begin()+1) {
|
---|
115 | histograms[6]->Fill(it->Pt());
|
---|
116 | histograms[7]->Fill(it->Eta());
|
---|
117 | histograms[8]->Fill(it->Phi());
|
---|
118 | }
|
---|
119 | }
|
---|
120 | }
|
---|
121 | Progress(int(100*i/entries));
|
---|
122 | }
|
---|
123 | treeReader->ReadEntry(event_id);
|
---|
124 | Progress(100);
|
---|
125 | }
|
---|
126 |
|
---|
127 | void DelphesPlotSummary::Draw() {
|
---|
128 | for(std::map< TString, TCanvas* >::iterator it=canvases_.begin(); it!=canvases_.end(); ++it) {
|
---|
129 | TCanvas* c = it->second;
|
---|
130 | std::vector<TH1F*> histograms = histograms_[it->first];
|
---|
131 | std::vector<TH1F*> eventProfiles = eventProfiles_[it->first];
|
---|
132 | std::vector<TMarker*> eventMarkers = eventMarkers_[it->first];
|
---|
133 | for(Int_t i=0;i<9;++i) {
|
---|
134 | c->cd(i+1);
|
---|
135 | if(histograms[i]->GetEntries()==0) continue;
|
---|
136 | histograms[i]->Draw();
|
---|
137 | if(i<3) {
|
---|
138 | eventProfiles[i]->Draw("same");
|
---|
139 | } else {
|
---|
140 | eventMarkers[i-3]->Draw("same");
|
---|
141 | }
|
---|
142 | }
|
---|
143 | c->Update();
|
---|
144 | }
|
---|
145 | }
|
---|
146 |
|
---|
147 | void DelphesPlotSummary::FillEvent() {
|
---|
148 | // clear event histograms and markers
|
---|
149 | for(std::map< TString, std::vector<TH1F*> >::iterator hv = eventProfiles_.begin(); hv!=eventProfiles_.end();++hv) {
|
---|
150 | for(std::vector<TH1F*>::iterator h = hv->second.begin(); h<hv->second.end();++h) {
|
---|
151 | (*h)->Reset();
|
---|
152 | }
|
---|
153 | }
|
---|
154 | for(std::map< TString, std::vector<TMarker*> >::iterator mv = eventMarkers_.begin(); mv!=eventMarkers_.end();++mv) {
|
---|
155 | for(std::vector<TMarker*>::iterator m = mv->second.begin(); m<mv->second.end();++m) {
|
---|
156 | (*m)->SetMarkerSize(0);
|
---|
157 | }
|
---|
158 | }
|
---|
159 | // loop over the elements and fill markers with event data
|
---|
160 | for(std::vector<DelphesBranchBase*>::iterator element = elements_->begin();element<elements_->end();++element) {
|
---|
161 | std::vector<TLorentzVector> vectors = (*element)->GetVectors();
|
---|
162 | std::sort(vectors.begin(), vectors.end(), vecsorter);
|
---|
163 | std::vector<TH1F*> hv = eventProfiles_[(*element)->GetName()];
|
---|
164 | TH1F* h1 = hv[0]; h1->Reset();
|
---|
165 | TH1F* h2 = hv[1]; h1->Reset();
|
---|
166 | TH1F* h3 = hv[2]; h1->Reset();
|
---|
167 | std::vector<TMarker*> mv = eventMarkers_[(*element)->GetName()];
|
---|
168 | for(std::vector<TLorentzVector>::iterator it=vectors.begin(); it<vectors.end();++it) {
|
---|
169 | h1->Fill(it->Pt());
|
---|
170 | h2->Fill(it->Eta());
|
---|
171 | h3->Fill(it->Phi());
|
---|
172 | if(it==vectors.begin()) {
|
---|
173 | mv[0]->SetX(it->Pt()); mv[0]->SetMarkerSize(3);
|
---|
174 | mv[1]->SetX(it->Eta()); mv[1]->SetMarkerSize(3);
|
---|
175 | mv[2]->SetX(it->Phi()); mv[2]->SetMarkerSize(3);
|
---|
176 | }
|
---|
177 | if(it==vectors.begin()+1) {
|
---|
178 | mv[3]->SetX(it->Pt()); mv[3]->SetMarkerSize(3);
|
---|
179 | mv[4]->SetX(it->Eta()); mv[4]->SetMarkerSize(3);
|
---|
180 | mv[5]->SetX(it->Phi()); mv[5]->SetMarkerSize(3);
|
---|
181 | }
|
---|
182 | }
|
---|
183 | }
|
---|
184 | }
|
---|