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) :
|
---|
26 | tab_(tab) {}
|
---|
27 |
|
---|
28 | DelphesPlotSummary::~DelphesPlotSummary() {}
|
---|
29 |
|
---|
30 | void DelphesPlotSummary::Progress(Int_t p)
|
---|
31 | {
|
---|
32 | Emit("Progress(Int_t)", p);
|
---|
33 | }
|
---|
34 |
|
---|
35 | void DelphesPlotSummary::Init(std::vector<DelphesBranchBase *> &elements)
|
---|
36 | {
|
---|
37 | elements_ = &elements;
|
---|
38 | // loop on the elements, and create tabs
|
---|
39 | for(std::vector<DelphesBranchBase *>::iterator data = elements.begin(); data < elements.end(); ++data)
|
---|
40 | {
|
---|
41 | // the canvas
|
---|
42 | TEveWindowSlot *slot = tab_->NewSlot();
|
---|
43 | TRootEmbeddedCanvas *trec = new TRootEmbeddedCanvas();
|
---|
44 | TCanvas *canvas = trec->GetCanvas();
|
---|
45 | TEveWindowFrame *wf = slot->MakeFrame(trec);
|
---|
46 | wf->SetElementName((*data)->GetName());
|
---|
47 | canvas->Divide(3, 3);
|
---|
48 | canvases_[(*data)->GetName()] = canvas;
|
---|
49 | // the histograms
|
---|
50 | TH1F *h;
|
---|
51 | std::vector<TH1F *> histograms;
|
---|
52 | histograms.reserve(9);
|
---|
53 | h = new TH1F(Form("%sPt", (*data)->GetName()), Form("%s Pt", (*data)->GetName()), 100, 0, -1);
|
---|
54 | histograms.push_back(h);
|
---|
55 | h = new TH1F(Form("%sEta", (*data)->GetName()), Form("%s Eta", (*data)->GetName()), 100, 0, -1);
|
---|
56 | histograms.push_back(h);
|
---|
57 | h = new TH1F(Form("%sPhi", (*data)->GetName()), Form("%s Phi", (*data)->GetName()), 100, 0, -1);
|
---|
58 | histograms.push_back(h);
|
---|
59 | h = new TH1F(Form("l%sPt", (*data)->GetName()), Form("leading %s Pt", (*data)->GetName()), 100, 0, -1);
|
---|
60 | histograms.push_back(h);
|
---|
61 | h = new TH1F(Form("l%sEta", (*data)->GetName()), Form("leading %s Eta", (*data)->GetName()), 100, 0, -1);
|
---|
62 | histograms.push_back(h);
|
---|
63 | h = new TH1F(Form("l%sPhi", (*data)->GetName()), Form("leading %s Phi", (*data)->GetName()), 100, 0, -1);
|
---|
64 | histograms.push_back(h);
|
---|
65 | h = new TH1F(Form("sl%sPt", (*data)->GetName()), Form("subleading %s Pt", (*data)->GetName()), 100, 0, -1);
|
---|
66 | histograms.push_back(h);
|
---|
67 | h = new TH1F(Form("sl%sEta", (*data)->GetName()), Form("subleading %s Eta", (*data)->GetName()), 100, 0, -1);
|
---|
68 | histograms.push_back(h);
|
---|
69 | h = new TH1F(Form("sl%sPhi", (*data)->GetName()), Form("subleading %s Phi", (*data)->GetName()), 100, 0, -1);
|
---|
70 | histograms.push_back(h);
|
---|
71 | histograms_[(*data)->GetName()] = histograms;
|
---|
72 | // the event histograms
|
---|
73 | TH1F *h1 = (TH1F *)histograms[0]->Clone();
|
---|
74 | h1->Reset();
|
---|
75 | h1->SetLineColor(kBlue);
|
---|
76 | TH1F *h2 = (TH1F *)histograms[1]->Clone();
|
---|
77 | h2->Reset();
|
---|
78 | h2->SetLineColor(kBlue);
|
---|
79 | TH1F *h3 = (TH1F *)histograms[2]->Clone();
|
---|
80 | h3->Reset();
|
---|
81 | h3->SetLineColor(kBlue);
|
---|
82 | std::vector<TH1F *> hv;
|
---|
83 | hv.push_back(h1);
|
---|
84 | hv.push_back(h2);
|
---|
85 | hv.push_back(h3);
|
---|
86 | eventProfiles_[(*data)->GetName()] = hv;
|
---|
87 | // the event markers
|
---|
88 | TMarker *m;
|
---|
89 | std::vector<TMarker *> mv;
|
---|
90 | m = new TMarker(0, 0, 29);
|
---|
91 | m->SetMarkerColor(kBlue);
|
---|
92 | m->SetMarkerSize(3);
|
---|
93 | mv.push_back(m);
|
---|
94 | m = new TMarker(0, 0, 29);
|
---|
95 | m->SetMarkerColor(kBlue);
|
---|
96 | m->SetMarkerSize(3);
|
---|
97 | mv.push_back(m);
|
---|
98 | m = new TMarker(0, 0, 29);
|
---|
99 | m->SetMarkerColor(kBlue);
|
---|
100 | m->SetMarkerSize(3);
|
---|
101 | mv.push_back(m);
|
---|
102 | m = new TMarker(0, 0, 29);
|
---|
103 | m->SetMarkerColor(kBlue);
|
---|
104 | m->SetMarkerSize(3);
|
---|
105 | mv.push_back(m);
|
---|
106 | m = new TMarker(0, 0, 29);
|
---|
107 | m->SetMarkerColor(kBlue);
|
---|
108 | m->SetMarkerSize(3);
|
---|
109 | mv.push_back(m);
|
---|
110 | m = new TMarker(0, 0, 29);
|
---|
111 | m->SetMarkerColor(kBlue);
|
---|
112 | m->SetMarkerSize(3);
|
---|
113 | mv.push_back(m);
|
---|
114 | eventMarkers_[(*data)->GetName()] = mv;
|
---|
115 | }
|
---|
116 | }
|
---|
117 |
|
---|
118 | void DelphesPlotSummary::FillSample(ExRootTreeReader *treeReader, Int_t event_id)
|
---|
119 | {
|
---|
120 | Int_t entries = treeReader->GetEntries();
|
---|
121 | for(Int_t i = 0; i < entries; ++i)
|
---|
122 | {
|
---|
123 | treeReader->ReadEntry(i);
|
---|
124 | for(std::vector<DelphesBranchBase *>::iterator element = elements_->begin(); element < elements_->end(); ++element)
|
---|
125 | {
|
---|
126 | std::vector<TLorentzVector> vectors = (*element)->GetVectors();
|
---|
127 | std::sort(vectors.begin(), vectors.end(), vecsorter);
|
---|
128 | std::vector<TH1F *> histograms = histograms_[(*element)->GetName()];
|
---|
129 | for(std::vector<TLorentzVector>::iterator it = vectors.begin(); it < vectors.end(); ++it)
|
---|
130 | {
|
---|
131 | histograms[0]->Fill(it->Pt());
|
---|
132 | histograms[1]->Fill(it->Eta());
|
---|
133 | histograms[2]->Fill(it->Phi());
|
---|
134 | if(it == vectors.begin())
|
---|
135 | {
|
---|
136 | histograms[3]->Fill(it->Pt());
|
---|
137 | histograms[4]->Fill(it->Eta());
|
---|
138 | histograms[5]->Fill(it->Phi());
|
---|
139 | }
|
---|
140 | if(it == vectors.begin() + 1)
|
---|
141 | {
|
---|
142 | histograms[6]->Fill(it->Pt());
|
---|
143 | histograms[7]->Fill(it->Eta());
|
---|
144 | histograms[8]->Fill(it->Phi());
|
---|
145 | }
|
---|
146 | }
|
---|
147 | }
|
---|
148 | Progress(int(100 * i / entries));
|
---|
149 | }
|
---|
150 | treeReader->ReadEntry(event_id);
|
---|
151 | Progress(100);
|
---|
152 | }
|
---|
153 |
|
---|
154 | void DelphesPlotSummary::Draw()
|
---|
155 | {
|
---|
156 | for(std::map<TString, TCanvas *>::iterator it = canvases_.begin(); it != canvases_.end(); ++it)
|
---|
157 | {
|
---|
158 | TCanvas *c = it->second;
|
---|
159 | std::vector<TH1F *> histograms = histograms_[it->first];
|
---|
160 | std::vector<TH1F *> eventProfiles = eventProfiles_[it->first];
|
---|
161 | std::vector<TMarker *> eventMarkers = eventMarkers_[it->first];
|
---|
162 | for(Int_t i = 0; i < 9; ++i)
|
---|
163 | {
|
---|
164 | c->cd(i + 1);
|
---|
165 | if(histograms[i]->GetEntries() == 0) continue;
|
---|
166 | histograms[i]->Draw();
|
---|
167 | if(i < 3)
|
---|
168 | {
|
---|
169 | eventProfiles[i]->Draw("same");
|
---|
170 | }
|
---|
171 | else
|
---|
172 | {
|
---|
173 | eventMarkers[i - 3]->Draw("same");
|
---|
174 | }
|
---|
175 | }
|
---|
176 | c->Update();
|
---|
177 | }
|
---|
178 | }
|
---|
179 |
|
---|
180 | void DelphesPlotSummary::FillEvent()
|
---|
181 | {
|
---|
182 | // clear event histograms and markers
|
---|
183 | for(std::map<TString, std::vector<TH1F *>>::iterator hv = eventProfiles_.begin(); hv != eventProfiles_.end(); ++hv)
|
---|
184 | {
|
---|
185 | for(std::vector<TH1F *>::iterator h = hv->second.begin(); h < hv->second.end(); ++h)
|
---|
186 | {
|
---|
187 | (*h)->Reset();
|
---|
188 | }
|
---|
189 | }
|
---|
190 | for(std::map<TString, std::vector<TMarker *>>::iterator mv = eventMarkers_.begin(); mv != eventMarkers_.end(); ++mv)
|
---|
191 | {
|
---|
192 | for(std::vector<TMarker *>::iterator m = mv->second.begin(); m < mv->second.end(); ++m)
|
---|
193 | {
|
---|
194 | (*m)->SetMarkerSize(0);
|
---|
195 | }
|
---|
196 | }
|
---|
197 | // loop over the elements and fill markers with event data
|
---|
198 | for(std::vector<DelphesBranchBase *>::iterator element = elements_->begin(); element < elements_->end(); ++element)
|
---|
199 | {
|
---|
200 | std::vector<TLorentzVector> vectors = (*element)->GetVectors();
|
---|
201 | std::sort(vectors.begin(), vectors.end(), vecsorter);
|
---|
202 | std::vector<TH1F *> hv = eventProfiles_[(*element)->GetName()];
|
---|
203 | TH1F *h1 = hv[0];
|
---|
204 | h1->Reset();
|
---|
205 | TH1F *h2 = hv[1];
|
---|
206 | h1->Reset();
|
---|
207 | TH1F *h3 = hv[2];
|
---|
208 | h1->Reset();
|
---|
209 | std::vector<TMarker *> mv = eventMarkers_[(*element)->GetName()];
|
---|
210 | for(std::vector<TLorentzVector>::iterator it = vectors.begin(); it < vectors.end(); ++it)
|
---|
211 | {
|
---|
212 | h1->Fill(it->Pt());
|
---|
213 | h2->Fill(it->Eta());
|
---|
214 | h3->Fill(it->Phi());
|
---|
215 | if(it == vectors.begin())
|
---|
216 | {
|
---|
217 | mv[0]->SetX(it->Pt());
|
---|
218 | mv[0]->SetMarkerSize(3);
|
---|
219 | mv[1]->SetX(it->Eta());
|
---|
220 | mv[1]->SetMarkerSize(3);
|
---|
221 | mv[2]->SetX(it->Phi());
|
---|
222 | mv[2]->SetMarkerSize(3);
|
---|
223 | }
|
---|
224 | if(it == vectors.begin() + 1)
|
---|
225 | {
|
---|
226 | mv[3]->SetX(it->Pt());
|
---|
227 | mv[3]->SetMarkerSize(3);
|
---|
228 | mv[4]->SetX(it->Eta());
|
---|
229 | mv[4]->SetMarkerSize(3);
|
---|
230 | mv[5]->SetX(it->Phi());
|
---|
231 | mv[5]->SetMarkerSize(3);
|
---|
232 | }
|
---|
233 | }
|
---|
234 | }
|
---|
235 | }
|
---|