Fork me on GitHub

source: git/examples/Validation.cpp@ 14fe596

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 14fe596 was ebf057e, checked in by Alexandre Mertens <alexandre.mertens@…>, 8 years ago

Merge branch 'dev_01' of https://github.com/delphes/delphes into ValidationWorkWithPavel

  • Property mode set to 100644
File size: 39.1 KB
Line 
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
20#include <iostream>
21#include <utility>
22#include <vector>
23#include <typeinfo>
24
25#include "TROOT.h"
26#include "TSystem.h"
27#include "TApplication.h"
28
29#include "TString.h"
30
31#include "TH1.h"
32#include "TH2.h"
33#include "TMath.h"
34#include "TStyle.h"
35#include "TGraph.h"
36#include "TCanvas.h"
37#include "THStack.h"
38#include "TLegend.h"
39#include "TPaveText.h"
40#include "TClonesArray.h"
41#include "TLorentzVector.h"
42#include "TGraphErrors.h"
43#include "TMultiGraph.h"
44
45#include "classes/DelphesClasses.h"
46
47#include "ExRootAnalysis/ExRootTreeReader.h"
48#include "ExRootAnalysis/ExRootTreeWriter.h"
49#include "ExRootAnalysis/ExRootTreeBranch.h"
50#include "ExRootAnalysis/ExRootResult.h"
51#include "ExRootAnalysis/ExRootUtilities.h"
52
53using namespace std;
54
55//------------------------------------------------------------------------------
56
57double ptrangemin = 10;
58double ptrangemax = 10000;
59static const int Nbins = 20;
60
61int objStyle = 1;
62int trackStyle = 7;
63int towerStyle = 3;
64
65Color_t objColor = kBlack;
66Color_t trackColor = kBlack;
67Color_t towerColor = kBlack;
68
69double effLegXmin = 0.22;
70double effLegXmax = 0.7;
71double effLegYmin = 0.22;
72double effLegYmax = 0.5;
73
74double resLegXmin = 0.62;
75double resLegXmax = 0.9;
76double resLegYmin = 0.52;
77double resLegYmax = 0.85;
78
79double topLeftLegXmin = 0.22;
80double topLeftLegXmax = 0.7;
81double topLeftLegYmin = 0.52;
82double topLeftLegYmax = 0.85;
83
84
85struct resolPlot
86{
87 TH1 *cenResolHist;
88 TH1 *fwdResolHist;
89 double ptmin;
90 double ptmax;
91 double xmin;
92 double xmax;
93 TString obj;
94
95 resolPlot();
96 resolPlot(double ptdown, double ptup, TString object);
97 void set(double ptdown, double ptup, TString object, double xmin = 0, double xmax = 2);
98 void print(){std::cout << ptmin << std::endl;}
99};
100
101
102resolPlot::resolPlot()
103{
104}
105
106resolPlot::resolPlot(double ptdown, double ptup, TString object)
107{
108 this->set(ptdown,ptup,object);
109}
110
111void resolPlot::set(double ptdown, double ptup, TString object, double xmin, double xmax)
112{
113 ptmin = ptdown;
114 ptmax = ptup;
115 obj = object;
116
117 cenResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", 200, xmin, xmax);
118 fwdResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", 200, xmin, xmax);
119}
120
121void HistogramsCollection(std::vector<resolPlot> *histos, double ptmin, double ptmax, TString obj, double xmin = 0, double xmax = 2)
122{
123 double width;
124 double ptdown;
125 double ptup;
126 resolPlot ptemp;
127
128 for (int i = 0; i < Nbins; i++)
129 {
130 width = (ptmax - ptmin) / Nbins;
131 ptdown = TMath::Power(10,ptmin + i * width );
132 ptup = TMath::Power(10,ptmin + (i+1) * width );
133 ptemp.set(ptdown, ptup, obj, xmin, xmax);
134 histos->push_back(ptemp);
135 }
136}
137
138//------------------------------------------------------------------------------
139
140class ExRootResult;
141class ExRootTreeReader;
142
143//------------------------------------------------------------------------------
144
145void BinLogX(TH1*h)
146{
147 TAxis *axis = h->GetXaxis();
148 int bins = axis->GetNbins();
149
150 Axis_t from = axis->GetXmin();
151 Axis_t to = axis->GetXmax();
152 Axis_t width = (to - from) / bins;
153 Axis_t *new_bins = new Axis_t[bins + 1];
154
155 for (int i = 0; i <= bins; i++)
156 {
157 new_bins[i] = TMath::Power(10, from + i * width);
158 }
159 axis->Set(bins, new_bins);
160 delete new_bins;
161}
162
163
164//------------------------------------------------------------------------------
165
166template<typename T>
167std::pair<TH1D*, TH1D*> GetEff(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, ExRootTreeReader *treeReader)
168{
169
170 cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
171
172 Long64_t allEntries = treeReader->GetEntries();
173
174 GenParticle *particle;
175 T *recoObj;
176
177 TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
178
179 Float_t deltaR;
180 Float_t pt, eta;
181 Long64_t entry;
182
183 Int_t i, j;
184
185 TH1D *histGenPtcen = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
186 TH1D *histRecoPtcen = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
187 TH1D *histGenPtfwd = new TH1D(name+" gen spectra Eta",name+" gen spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
188 TH1D *histRecoPtfwd = new TH1D(name+" reco spectra Eta",name+" reco spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
189
190 histGenPtcen->SetDirectory(0);
191 histRecoPtcen->SetDirectory(0);
192 histGenPtfwd->SetDirectory(0);
193 histRecoPtfwd->SetDirectory(0);
194
195 BinLogX(histGenPtcen);
196 BinLogX(histRecoPtcen);
197 BinLogX(histGenPtfwd);
198 BinLogX(histRecoPtfwd);
199
200 // Loop over all events
201 for(entry = 0; entry < allEntries; ++entry)
202 {
203 // Load selected branches with data from specified event
204 treeReader->ReadEntry(entry);
205
206 // Loop over all generated particle in event
207 for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
208 {
209
210 particle = (GenParticle*) branchParticle->At(i);
211 genMomentum = particle->P4();
212
213 deltaR = 999;
214
215 if (particle->PID == pdgID && genMomentum.Pt() > ptrangemin && genMomentum.Pt() < ptrangemax )
216 {
217
218 // Loop over all reco object in event
219 for(j = 0; j < branchReco->GetEntriesFast(); ++j)
220 {
221 recoObj = (T*)branchReco->At(j);
222 recoMomentum = recoObj->P4();
223 // this is simply to avoid warnings from initial state particle
224 // having infite rapidity ...
225 //if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
226
227 // take the closest parton candidate
228 if(TMath::Abs(pdgID) == 5)
229 {
230 Jet *jet = (Jet *)recoObj;
231 if(jet->BTag != 1) continue;
232 }
233 if(TMath::Abs(pdgID) == 15)
234 {
235 Jet *jet = (Jet *)recoObj;
236 if(jet->TauTag != 1) continue;
237 }
238 if(genMomentum.DeltaR(recoMomentum) < deltaR)
239 {
240 deltaR = genMomentum.DeltaR(recoMomentum);
241 bestRecoMomentum = recoMomentum;
242 }
243 }
244
245 pt = genMomentum.Pt();
246 eta = genMomentum.Eta();
247
248 if (TMath::Abs(eta) < 1.5)
249 {
250 histGenPtcen->Fill(pt);
251 if(deltaR < 0.3) { histRecoPtcen->Fill(pt); }
252 }
253 else if (TMath::Abs(eta) < 2.5)
254 {
255 histGenPtfwd->Fill(pt);
256 if(deltaR < 0.3) { histRecoPtfwd->Fill(pt); }
257
258 }
259 }
260 }
261 }
262
263
264 std::pair<TH1D*,TH1D*> histos;
265
266 histRecoPtcen->Divide(histGenPtcen);
267 histRecoPtfwd->Divide(histGenPtfwd);
268
269 histos.first = histRecoPtcen;
270 histos.second = histRecoPtfwd;
271
272 return histos;
273}
274
275template<typename T>
276void GetEres(std::vector<resolPlot> *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, ExRootTreeReader *treeReader)
277{
278 Long64_t allEntries = treeReader->GetEntries();
279
280 cout << "** Computing resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
281
282 GenParticle *particle;
283 T* recoObj;
284
285 TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
286
287 Float_t deltaR;
288 Float_t pt, eta;
289 Long64_t entry;
290
291 Int_t i, j, bin;
292
293 // Loop over all events
294 for(entry = 0; entry < allEntries; ++entry)
295 {
296 // Load selected branches with data from specified event
297 treeReader->ReadEntry(entry);
298
299 // Loop over all reconstructed jets in event
300 for(i = 0; i < branchReco->GetEntriesFast(); ++i)
301 {
302 recoObj = (T*) branchReco->At(i);
303 recoMomentum = recoObj->P4();
304
305 deltaR = 999;
306
307 // Loop over all hard partons in event
308 for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
309 {
310 particle = (GenParticle*) branchParticle->At(j);
311 if (particle->PID == pdgID && particle->Status == 1)
312 {
313 genMomentum = particle->P4();
314
315 // this is simply to avoid warnings from initial state particle
316 // having infite rapidity ...
317 if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
318
319 // take the closest parton candidate
320 if(genMomentum.DeltaR(recoMomentum) < deltaR)
321 {
322 deltaR = genMomentum.DeltaR(recoMomentum);
323 bestGenMomentum = genMomentum;
324 }
325 }
326 }
327
328 if(deltaR < 0.3)
329 {
330 pt = bestGenMomentum.E();
331 eta = TMath::Abs(bestGenMomentum.Eta());
332
333 for (bin = 0; bin < Nbins; bin++)
334 {
335 if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 2.5)
336 {
337 if (eta < 1.5) {histos->at(bin).cenResolHist->Fill(recoMomentum.E()/bestGenMomentum.E());}
338 else if (eta < 2.5) {histos->at(bin).fwdResolHist->Fill(recoMomentum.E()/bestGenMomentum.E());}
339 }
340 }
341 }
342 }
343 }
344}
345void GetJetsEres(std::vector<resolPlot> *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader)
346{
347
348 Long64_t allEntries = treeReader->GetEntries();
349
350 cout << "** Computing resolution of " << branchJet->GetName() << " induced by " << branchGenJet->GetName() << endl;
351
352 Jet *jet, *genjet;
353
354 TLorentzVector jetMomentum, genJetMomentum, bestGenJetMomentum;
355
356 Float_t deltaR;
357 Float_t pt, eta;
358 Long64_t entry;
359
360 Int_t i, j, bin;
361
362 // Loop over all events
363 for(entry = 0; entry < allEntries; ++entry)
364 {
365 // Load selected branches with data from specified event
366 treeReader->ReadEntry(entry);
367
368 if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
369
370 // Loop over all reconstructed jets in event
371 for(i = 0; i < TMath::Min(2,branchJet->GetEntriesFast()); ++i) //branchJet->GetEntriesFast(); ++i)
372 {
373
374 jet = (Jet*) branchJet->At(i);
375 jetMomentum = jet->P4();
376
377 deltaR = 999;
378
379 // Loop over all hard partons in event
380 for(j = 0; j < TMath::Min(2,branchGenJet->GetEntriesFast()); ++j)
381 {
382 genjet = (Jet*) branchGenJet->At(j);
383
384 genJetMomentum = genjet->P4();
385
386 // this is simply to avoid warnings from initial state particle
387 // having infite rapidity ...
388 if(genJetMomentum.Px() == 0 && genJetMomentum.Py() == 0) continue;
389
390 // take the closest parton candidate
391 if(genJetMomentum.DeltaR(jetMomentum) < deltaR)
392 {
393 deltaR = genJetMomentum.DeltaR(jetMomentum);
394 bestGenJetMomentum = genJetMomentum;
395 }
396 }
397
398 if(deltaR < 0.25)
399 {
400 pt = genJetMomentum.Pt();
401 eta = TMath::Abs(genJetMomentum.Eta());
402
403 for (bin = 0; bin < Nbins; bin++)
404 {
405 if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 1.5)
406 {
407 histos->at(bin).cenResolHist->Fill(jetMomentum.Pt()/bestGenJetMomentum.Pt());
408 }
409 else if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 2.5)
410 {
411 histos->at(bin).fwdResolHist->Fill(jetMomentum.Pt()/bestGenJetMomentum.Pt());
412 }
413
414 }
415 }
416 }
417 }
418}
419
420void GetMetres(std::vector<resolPlot> *histos, TClonesArray *branchScalarHT, TClonesArray *branchMet, TClonesArray *branchJet, ExRootTreeReader *treeReader)
421{
422
423 Long64_t allEntries = treeReader->GetEntries();
424
425 cout << "** Computing resolution of " << branchMet->GetName() << " vs " << branchScalarHT->GetName() << endl;
426
427 MissingET *met;
428 ScalarHT *scalarHT;
429
430 Long64_t entry;
431
432 Int_t bin;
433 Double_t ht;
434
435 Jet *jet;
436 TLorentzVector p1, p2;
437
438 // Loop over all events
439 for(entry = 0; entry < allEntries; ++entry)
440 {
441 // Load selected branches with data from specified event
442 treeReader->ReadEntry(entry);
443
444 if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
445
446 if (branchJet->GetEntriesFast() > 1)
447 {
448
449 jet = (Jet*) branchJet->At(0);
450 p1 = jet->P4();
451 jet = (Jet*) branchJet->At(1);
452 p2 = jet->P4();
453
454 met = (MissingET*) branchMet->At(0);
455 scalarHT = (ScalarHT*) branchScalarHT->At(0);
456 ht = scalarHT->HT;
457
458 if(p1.Pt() < 0.75*ht/2) continue;
459 if(p2.Pt() < 0.75*ht/2) continue;
460
461 for (bin = 0; bin < Nbins; bin++)
462 {
463 if(ht > histos->at(bin).ptmin && ht < histos->at(bin).ptmax )
464 {
465 histos->at(bin).cenResolHist->Fill(met->P4().Px());
466 histos->at(bin).fwdResolHist->Fill(met->P4().Py());
467 }
468 }
469 }
470 }
471}
472
473
474std::pair<Double_t, Double_t> GausFit(TH1* hist)
475{
476 TF1 *f1 = new TF1("f1", "gaus", hist->GetMean()-2*hist->GetRMS(), hist->GetMean()+2*hist->GetRMS());
477 hist->Fit("f1","RQ");
478
479 TF1 *f2 = new TF1("f2", "gaus", f1->GetParameter(1) - 2*f1->GetParameter(2), f1->GetParameter(1) + 2*f1->GetParameter(2));
480 hist->Fit("f2","RQ");
481
482 Double_t sig = f2->GetParameter(2);
483 Double_t sigErr = f2->GetParError(2);
484
485 delete f1;
486 return make_pair (sig, sigErr);
487}
488
489
490TGraphErrors EresGraph(std::vector<resolPlot> *histos, bool central, bool rms = false)
491{
492 Int_t bin;
493 Int_t count = 0;
494 TGraphErrors gr = TGraphErrors(Nbins/2);
495 Double_t sig = 0;
496 Double_t sigErr = 0;
497 for (bin = 0; bin < Nbins; bin++)
498 {
499 if (central == true && histos->at(bin).cenResolHist->GetEntries() > 100)
500 {
501 std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).cenResolHist);
502 if (rms == true)
503 {
504 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
505 gr.SetPointError(count,0, sigvalues.second); // to correct
506 }
507 else
508 {
509 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
510 gr.SetPointError(count,0, sigvalues.second);
511 }
512 count++;
513 }
514
515 else if (central == false && histos->at(bin).fwdResolHist->GetEntries() > 10)
516 {
517 std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).fwdResolHist);
518 if (rms == true)
519 {
520 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
521 gr.SetPointError(count,0, sigvalues.second); // to correct
522 }
523 else
524 {
525 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
526 gr.SetPointError(count,0, sigvalues.second);
527 }
528 count++;
529 }
530
531 }
532 return gr;
533}
534
535
536//------------------------------------------------------------------------------
537
538
539// type 1 : object, 2 : track, 3 : tower
540
541void addGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int type)
542{
543
544 gr->SetLineWidth(2);
545
546 switch ( type )
547 {
548 case 1:
549 gr->SetLineColor(objColor);
550 gr->SetLineStyle(objStyle);
551 std::cout << "Adding " << gr->GetName() << std::endl;
552 mg->Add(gr);
553 leg->AddEntry(gr,"Reco","l");
554 break;
555
556 case 2:
557 gr->SetLineColor(trackColor);
558 gr->SetLineStyle(trackStyle);
559 mg->Add(gr);
560 leg->AddEntry(gr,"Track","l");
561 break;
562
563 case 3:
564 gr->SetLineColor(towerColor);
565 gr->SetLineStyle(towerStyle);
566 mg->Add(gr);
567 leg->AddEntry(gr,"Tower","l");
568 break;
569
570 case 0:
571 gr->SetLineColor(objColor);
572 gr->SetLineStyle(objStyle);
573 mg->Add(gr);
574 break;
575
576 default:
577 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
578 break;
579 }
580}
581
582void addHist(TH1D *h, TLegend *leg, int type)
583{
584 h->SetLineWidth(2);
585
586 switch ( type )
587 {
588 case 1:
589 h->SetLineColor(objColor);
590 h->SetLineStyle(objStyle);
591 leg->AddEntry(h,"Reco","l");
592 break;
593
594 case 2:
595 h->SetLineColor(trackColor);
596 h->SetLineStyle(trackStyle);
597 leg->AddEntry(h,"Track","l");
598 break;
599
600 case 3:
601 h->SetLineColor(towerColor);
602 h->SetLineStyle(towerStyle);
603 leg->AddEntry(h,"Tower","l");
604 break;
605
606 case 0:
607 h->SetLineColor(objColor);
608 h->SetLineStyle(objStyle);
609 break;
610
611 default:
612 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
613 break;
614 }
615}
616
617void DrawAxis(TMultiGraph *mg, TLegend *leg, double max)
618{
619 mg->SetMinimum(0.);
620 mg->SetMaximum(max);
621 mg->GetXaxis()->SetLimits(ptrangemin,ptrangemax);
622 mg->GetYaxis()->SetTitle("resolution");
623 mg->GetXaxis()->SetTitle("p_{T} [GeV]");
624 mg->GetYaxis()->SetTitleSize(0.07);
625 mg->GetXaxis()->SetTitleSize(0.07);
626 mg->GetYaxis()->SetLabelSize(0.06);
627 mg->GetXaxis()->SetLabelSize(0.06);
628 mg->GetYaxis()->SetLabelOffset(0.03);
629 mg->GetYaxis()->SetTitleOffset(1.4);
630 mg->GetXaxis()->SetTitleOffset(1.4);
631
632 mg->GetYaxis()->SetNdivisions(505);
633
634 leg->SetBorderSize(0);
635 leg->SetShadowColor(0);
636 leg->SetFillColor(0);
637 leg->SetFillStyle(0);
638
639 gStyle->SetOptTitle(0);
640 gPad->SetLogx();
641 gPad->SetBottomMargin(0.2);
642 gPad->SetLeftMargin(0.2);
643 gPad->Modified();
644 gPad->Update();
645
646}
647
648void DrawAxis(TH1D *h, TLegend *leg, int type)
649{
650
651 h->GetYaxis()->SetRangeUser(0,1.0);
652 if (type == 0) h->GetXaxis()->SetTitle("p_{T} [GeV]");
653 else h->GetXaxis()->SetTitle("#eta");
654 h->GetYaxis()->SetTitle("efficiency");
655 h->GetYaxis()->SetTitleSize(0.07);
656 h->GetXaxis()->SetTitleSize(0.07);
657 h->GetYaxis()->SetLabelSize(0.06);
658 h->GetXaxis()->SetLabelSize(0.06);
659 h->GetYaxis()->SetLabelOffset(0.03);
660 h->GetYaxis()->SetTitleOffset(1.3);
661 h->GetXaxis()->SetTitleOffset(1.4);
662
663 h->GetYaxis()->SetNdivisions(505);
664
665 leg->SetBorderSize(0);
666 leg->SetShadowColor(0);
667 leg->SetFillColor(0);
668 leg->SetFillStyle(0);
669
670 gStyle->SetOptTitle(0);
671 gStyle->SetOptStat(0);
672 gPad->SetBottomMargin(0.2);
673 gPad->SetLeftMargin(0.2);
674
675 gPad->Modified();
676 gPad->Update();
677
678}
679
680
681void Validation(const char *inputFileElectron, const char *inputFileMuon, const char *inputFilePhoton, const char *inputFileJet, const char *inputFileBJet, const char *inputFileTauJet, const char *outputFile)
682{
683 TChain *chainElectron = new TChain("Delphes");
684 chainElectron->Add(inputFileElectron);
685 ExRootTreeReader *treeReaderElectron = new ExRootTreeReader(chainElectron);
686
687 TChain *chainMuon = new TChain("Delphes");
688 chainMuon->Add(inputFileMuon);
689 ExRootTreeReader *treeReaderMuon = new ExRootTreeReader(chainMuon);
690
691 TChain *chainPhoton = new TChain("Delphes");
692 chainPhoton->Add(inputFilePhoton);
693 ExRootTreeReader *treeReaderPhoton = new ExRootTreeReader(chainPhoton);
694
695 TChain *chainJet = new TChain("Delphes");
696 chainJet->Add(inputFileJet);
697 ExRootTreeReader *treeReaderJet = new ExRootTreeReader(chainJet);
698
699 TChain *chainBJet = new TChain("Delphes");
700 chainBJet->Add(inputFileBJet);
701 ExRootTreeReader *treeReaderBJet = new ExRootTreeReader(chainBJet);
702
703 TChain *chainTauJet = new TChain("Delphes");
704 chainTauJet->Add(inputFileTauJet);
705 ExRootTreeReader *treeReaderTauJet = new ExRootTreeReader(chainTauJet);
706
707 TClonesArray *branchParticleElectron = treeReaderElectron->UseBranch("Particle");
708 TClonesArray *branchTrackElectron = treeReaderElectron->UseBranch("Track");
709 TClonesArray *branchTowerElectron = treeReaderElectron->UseBranch("Tower");
710 TClonesArray *branchElectron = treeReaderElectron->UseBranch("Electron");
711
712 TClonesArray *branchParticleMuon = treeReaderMuon->UseBranch("Particle");
713 TClonesArray *branchTrackMuon = treeReaderMuon->UseBranch("Track");
714 TClonesArray *branchMuon = treeReaderMuon->UseBranch("Muon");
715
716 TClonesArray *branchParticlePhoton = treeReaderPhoton->UseBranch("Particle");
717 TClonesArray *branchTowerPhoton = treeReaderPhoton->UseBranch("Tower");
718 TClonesArray *branchPhoton = treeReaderPhoton->UseBranch("Photon");
719
720 TClonesArray *branchGenJet = treeReaderJet->UseBranch("GenJet");
721 TClonesArray *branchPFJet = treeReaderJet->UseBranch("Jet");
722 TClonesArray *branchCaloJet = treeReaderJet->UseBranch("CaloJet");
723
724 TClonesArray *branchParticleBJet = treeReaderBJet->UseBranch("Particle");
725 TClonesArray *branchPFBJet = treeReaderBJet->UseBranch("Jet");
726
727 TClonesArray *branchParticleTauJet = treeReaderTauJet->UseBranch("Particle");
728 TClonesArray *branchPFTauJet = treeReaderTauJet->UseBranch("Jet");
729
730 TClonesArray *branchScalarHT = treeReaderJet->UseBranch("ScalarHT");
731 TClonesArray *branchMet = treeReaderJet->UseBranch("MissingET");
732
733 ///////////////
734 // Electrons //
735 ///////////////
736
737 // Reconstruction efficiency
738 TString elecs = "Electron";
739 int elID = 11;
740 std::pair<TH1D*,TH1D*> histos_el = GetEff<Electron>(branchElectron, branchParticleElectron, "Electron", elID, treeReaderElectron);
741
742 // tracking reconstruction efficiency
743 std::pair <TH1D*,TH1D*> histos_eltrack = GetEff<Track>(branchTrackElectron, branchParticleElectron, "electronTrack", elID, treeReaderElectron);
744
745 // Tower reconstruction efficiency
746 std::pair <TH1D*,TH1D*> histos_eltower = GetEff<Tower>(branchTowerElectron, branchParticleElectron, "electronTower", elID, treeReaderElectron);
747
748 // Electron Energy Resolution
749 std::vector<resolPlot> plots_el;
750 HistogramsCollection(&plots_el, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electrons");
751 GetEres<Electron>(&plots_el, branchElectron, branchParticleElectron, elID, treeReaderElectron);
752 TGraphErrors gr_el = EresGraph(&plots_el, true);
753 TGraphErrors gr_elFwd = EresGraph(&plots_el, false);
754 gr_el.SetName("Electron");
755 gr_elFwd.SetName("ElectronFwd");
756
757 // Electron Track Energy Resolution
758 std::vector<resolPlot> plots_eltrack;
759 HistogramsCollection(&plots_eltrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTracks");
760 GetEres<Track>(&plots_eltrack, branchTrackElectron, branchParticleElectron, elID, treeReaderElectron);
761 TGraphErrors gr_eltrack = EresGraph(&plots_eltrack, true);
762 TGraphErrors gr_eltrackFwd = EresGraph(&plots_eltrack, false);
763 gr_eltrack.SetName("ElectronTracks");
764 gr_eltrackFwd.SetName("ElectronTracksFwd");
765
766 // Electron Tower Energy Resolution
767 std::vector<resolPlot> plots_eltower;
768 HistogramsCollection(&plots_eltower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTower");
769 GetEres<Tower>(&plots_eltower, branchTowerElectron, branchParticleElectron, elID, treeReaderElectron);
770 TGraphErrors gr_eltower = EresGraph(&plots_eltower, true);
771 TGraphErrors gr_eltowerFwd = EresGraph(&plots_eltower, false);
772 gr_eltower.SetName("ElectronTower");
773 gr_eltrackFwd.SetName("ElectronTracksFwd");
774
775 // Canvases
776 TString elEff = "electronEff";
777 TCanvas *C_el1 = new TCanvas(elEff,elEff, 1600, 600);
778 C_el1->Divide(2);
779 C_el1->cd(1);
780 TLegend *leg_el1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
781 leg_el1->SetHeader("#splitline{electrons}{|#eta| < 1.5}");
782 leg_el1->AddEntry("","","");
783
784 gPad->SetLogx();
785 histos_eltrack.first->Draw("][");
786 addHist(histos_eltrack.first, leg_el1, 2);
787 histos_el.first->Draw("same ][");
788 addHist(histos_el.first, leg_el1, 1);
789 DrawAxis(histos_eltrack.first, leg_el1, 0);
790
791 leg_el1->Draw();
792
793 C_el1->cd(2);
794 TLegend *leg_el2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
795 leg_el2->SetHeader("#splitline{electrons}{1.5 < |#eta| < 2.5}");
796 leg_el2->AddEntry("","","");
797
798 gPad->SetLogx();
799 histos_eltrack.second->Draw("][");
800 addHist(histos_eltrack.second, leg_el2, 2);
801 histos_el.second->Draw("same ][");
802 addHist(histos_el.second, leg_el2, 1);
803
804 DrawAxis(histos_eltrack.second, leg_el2, 0);
805 leg_el2->Draw();
806
807 TString elRes = "electronERes";
808 TString elResFwd = "electronEResForward";
809 TCanvas *C_el2 = new TCanvas(elRes,elRes, 1600, 600);
810 C_el2->Divide(2);
811 C_el2->cd(1);
812 TMultiGraph *mg_el = new TMultiGraph(elRes,elRes);
813 TLegend *leg_el = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
814 leg_el->SetHeader("#splitline{electrons}{|#eta| < 1.5}");
815 leg_el->AddEntry("","","");
816
817 addGraph(mg_el, &gr_eltower, leg_el, 3);
818 addGraph(mg_el, &gr_eltrack, leg_el, 2);
819 addGraph(mg_el, &gr_el, leg_el, 1);
820
821 mg_el->Draw("ACX");
822 leg_el->Draw();
823
824 DrawAxis(mg_el, leg_el, 0.1);
825
826 C_el2->cd(2);
827 TMultiGraph *mg_elFwd = new TMultiGraph(elResFwd,elResFwd);
828 TLegend *leg_elFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
829 leg_elFwd->SetHeader("#splitline{electrons}{1.5 < |#eta| < 2.5}");
830 leg_elFwd->AddEntry("","","");
831
832 addGraph(mg_elFwd, &gr_eltowerFwd, leg_elFwd, 3);
833 addGraph(mg_elFwd, &gr_eltrackFwd, leg_elFwd, 2);
834 addGraph(mg_elFwd, &gr_elFwd, leg_elFwd, 1);
835
836 mg_elFwd->Draw("ACX");
837 leg_elFwd->Draw();
838
839 DrawAxis(mg_elFwd, leg_elFwd, 0.2);
840
841 C_el1->Print("delphes_validation.pdf(","pdf");
842 C_el2->Print("delphes_validation.pdf","pdf");
843
844 gDirectory->cd(0);
845
846 ///////////
847 // Muons //
848 ///////////
849
850 // Reconstruction efficiency
851 int muID = 13;
852 std::pair<TH1D*,TH1D*> histos_mu = GetEff<Muon>(branchMuon, branchParticleMuon,"Muon", muID, treeReaderMuon);
853
854 // muon tracking reconstruction efficiency
855 std::pair <TH1D*,TH1D*> histos_mutrack = GetEff<Track>(branchTrackMuon, branchParticleMuon, "muonTrack", muID, treeReaderMuon);
856
857 // Muon Energy Resolution
858 std::vector<resolPlot> plots_mu;
859 HistogramsCollection(&plots_mu, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muons");
860 GetEres<Muon>(&plots_mu, branchMuon, branchParticleMuon, muID, treeReaderMuon);
861 TGraphErrors gr_mu = EresGraph(&plots_mu, true);
862 TGraphErrors gr_muFwd = EresGraph(&plots_mu, false);
863 gr_mu.SetName("Muon");
864 gr_muFwd.SetName("MuonFwd");
865
866 // Muon Track Energy Resolution
867 std::vector<resolPlot> plots_mutrack;
868 HistogramsCollection(&plots_mutrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muonsTracks");
869 GetEres<Track>(&plots_mutrack, branchTrackMuon, branchParticleMuon, muID, treeReaderMuon);
870 TGraphErrors gr_mutrack = EresGraph(&plots_mutrack, true);
871 TGraphErrors gr_mutrackFwd = EresGraph(&plots_mutrack, false);
872 gr_mutrackFwd.SetName("MuonTracksFwd");
873
874 // Canvas
875
876 TString muEff = "muonEff";
877 TCanvas *C_mu1 = new TCanvas(muEff,muEff, 1600, 600);
878 C_mu1->Divide(2);
879 C_mu1->cd(1);
880 TLegend *leg_mu1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
881 leg_mu1->SetHeader("#splitline{muons}{|#eta| < 1.5}");
882 leg_mu1->AddEntry("","","");
883
884
885 gPad->SetLogx();
886 histos_mutrack.first->Draw("][");
887 addHist(histos_mutrack.first, leg_mu1, 2);
888 histos_mu.first->Draw("same ][");
889 addHist(histos_mu.first, leg_mu1, 1);
890
891 DrawAxis(histos_mutrack.first, leg_mu1, 0);
892
893 leg_mu1->Draw();
894
895 C_mu1->cd(2);
896 TLegend *leg_mu2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
897 leg_mu2->SetHeader("#splitline{muons}{1.5 < |#eta| < 2.5}");
898 leg_mu2->AddEntry("","","");
899
900 gPad->SetLogx();
901 histos_mutrack.second->Draw("][");
902 addHist(histos_mutrack.second, leg_mu2, 2);
903 histos_mu.second->Draw("same ][");
904 addHist(histos_mu.second, leg_mu2, 1);
905
906 DrawAxis(histos_mutrack.second, leg_mu2, 1);
907 leg_mu2->Draw();
908
909 TString muRes = "muonERes";
910 TString muResFwd = "muonEResFwd";
911
912 TCanvas *C_mu = new TCanvas(muRes,muRes, 1600, 600);
913 C_mu->Divide(2);
914 C_mu->cd(1);
915 TMultiGraph *mg_mu = new TMultiGraph(muRes,muRes);
916 TLegend *leg_mu = new TLegend(topLeftLegXmin,topLeftLegYmin,topLeftLegXmax,topLeftLegYmax);
917 leg_mu->SetHeader("#splitline{muons}{|#eta| < 1.5}");
918 leg_mu->AddEntry("","","");
919
920 addGraph(mg_mu, &gr_mutrack, leg_mu, 2);
921 addGraph(mg_mu, &gr_mu, leg_mu, 1);
922
923 mg_mu->Draw("ACX");
924 leg_mu->Draw();
925
926 DrawAxis(mg_mu, leg_mu, 0.3);
927
928 C_mu->cd(2);
929 TMultiGraph *mg_muFwd = new TMultiGraph(muResFwd,muResFwd);
930 TLegend *leg_muFwd = new TLegend(topLeftLegXmin,topLeftLegYmin,topLeftLegXmax,topLeftLegYmax);
931 leg_muFwd->SetHeader("#splitline{muons}{1.5 < |#eta| < 2.5}");
932 leg_muFwd->AddEntry("","","");
933
934 addGraph(mg_muFwd, &gr_mutrackFwd, leg_muFwd, 2);
935 addGraph(mg_muFwd, &gr_muFwd, leg_muFwd, 1);
936
937 mg_muFwd->Draw("ACX");
938 leg_muFwd->Draw();
939
940 DrawAxis(mg_muFwd, leg_muFwd, 0.3);
941
942 C_mu1->Print("delphes_validation.pdf","pdf");
943 C_mu->Print("delphes_validation.pdf","pdf");
944
945 gDirectory->cd(0);
946
947 /////////////
948 // Photons //
949 /////////////
950
951 // Reconstruction efficiency
952 int phID = 22;
953 std::pair<TH1D*,TH1D*> histos_ph = GetEff<Electron>(branchPhoton, branchParticlePhoton, "Photon", phID, treeReaderPhoton);
954 std::pair<TH1D*,TH1D*> histos_phtower = GetEff<Electron>(branchTowerPhoton, branchParticlePhoton, "Photon", phID, treeReaderPhoton);
955
956 // Photon Energy Resolution
957 std::vector<resolPlot> plots_ph;
958 HistogramsCollection(&plots_ph, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photons");
959 GetEres<Photon>(&plots_ph, branchPhoton, branchParticlePhoton, phID, treeReaderPhoton);
960 TGraphErrors gr_ph = EresGraph(&plots_ph, true);
961 TGraphErrors gr_phFwd = EresGraph(&plots_ph, false);
962 gr_ph.SetName("Photon");
963 gr_phFwd.SetName("PhotonFwd");
964
965
966 // Photon Tower Energy Resolution
967 std::vector<resolPlot> plots_phtower;
968 HistogramsCollection(&plots_phtower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photonsTower");
969 GetEres<Tower>(&plots_phtower, branchTowerPhoton, branchParticlePhoton, phID, treeReaderPhoton);
970 TGraphErrors gr_phtower = EresGraph(&plots_phtower, true);
971 TGraphErrors gr_phtowerFwd = EresGraph(&plots_phtower, false);
972 gr_phtower.SetName("PhotonTower");
973 gr_phtowerFwd.SetName("PhotonTowerFwd");
974
975 // Canvas
976
977 TString phEff = "photonEff";
978 TCanvas *C_ph1 = new TCanvas(phEff,phEff, 1600, 600);
979 C_ph1->Divide(2);
980 C_ph1->cd(1);
981 TLegend *leg_ph1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
982 leg_ph1->SetHeader("#splitline{photons}{|#eta| < 1.5}");
983 leg_ph1->AddEntry("","","");
984
985
986 gPad->SetLogx();
987 histos_phtower.first->Draw("][");
988 addHist(histos_phtower.first, leg_ph1, 3);
989 histos_ph.first->Draw("same ][");
990 addHist(histos_ph.first, leg_ph1, 1);
991
992 DrawAxis(histos_phtower.first, leg_ph1, 0);
993 leg_ph1->Draw();
994
995 C_ph1->cd(2);
996 TLegend *leg_ph2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
997 leg_ph2->SetHeader("#splitline{photons}{1.5 < |#eta| < 2.5}");
998 leg_ph2->AddEntry("","","");
999
1000
1001 gPad->SetLogx();
1002 histos_phtower.second->Draw("][");
1003 addHist(histos_phtower.second, leg_ph2, 3);
1004 histos_ph.second->Draw("same ][");
1005 addHist(histos_ph.second, leg_ph2, 1);
1006
1007 DrawAxis(histos_phtower.second, leg_ph2, 1);
1008 leg_ph2->Draw();
1009
1010 TString phRes = "phERes";
1011 TString phResFwd = "phEResFwd";
1012
1013 TCanvas *C_ph = new TCanvas(phRes,phRes, 1600, 600);
1014 C_ph->Divide(2);
1015 C_ph->cd(1);
1016 TMultiGraph *mg_ph = new TMultiGraph(phRes,phRes);
1017 TLegend *leg_ph = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1018 leg_ph->SetHeader("#splitline{photons}{|#eta| < 1.5}");
1019 leg_ph->AddEntry("","","");
1020
1021 addGraph(mg_ph, &gr_phtower, leg_ph, 3);
1022 addGraph(mg_ph, &gr_ph, leg_ph, 1);
1023
1024 mg_ph->Draw("ACX");
1025 leg_ph->Draw();
1026
1027 DrawAxis(mg_ph, leg_ph, 0.3);
1028
1029 C_ph->cd(2);
1030 TMultiGraph *mg_phFwd = new TMultiGraph(phResFwd,phResFwd);
1031 TLegend *leg_phFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1032 leg_phFwd->SetHeader("#splitline{photons}{1.5 < |#eta| < 2.5}");
1033 leg_phFwd->AddEntry("","","");
1034
1035 addGraph(mg_phFwd, &gr_phtowerFwd, leg_phFwd, 3);
1036 addGraph(mg_phFwd, &gr_phFwd, leg_phFwd, 1);
1037
1038 mg_phFwd->Draw("ACX");
1039 leg_phFwd->Draw();
1040
1041 DrawAxis(mg_phFwd, leg_phFwd, 0.3);
1042
1043 C_ph1->Print("delphes_validation.pdf","pdf");
1044 C_ph->Print("delphes_validation.pdf","pdf");
1045
1046 gDirectory->cd(0);
1047
1048 //////////
1049 // Jets //
1050 //////////
1051
1052 // BJets Reconstruction efficiency
1053 int bID = 5;
1054 std::pair<TH1D*,TH1D*> histos_btag = GetEff<Jet>(branchPFBJet, branchParticleBJet,"BTag", bID, treeReaderBJet);
1055
1056 // TauJets Reconstruction efficiency
1057 int tauID = 15;
1058 std::pair<TH1D*,TH1D*> histos_tautag = GetEff<Jet>(branchPFTauJet, branchParticleTauJet,"TauTag", tauID, treeReaderTauJet);
1059
1060 // PFJets Energy Resolution
1061 std::vector<resolPlot> plots_pfjets;
1062 HistogramsCollection(&plots_pfjets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "PFJet");
1063 GetJetsEres(&plots_pfjets, branchPFJet, branchGenJet, treeReaderJet);
1064 TGraphErrors gr_pfjets = EresGraph(&plots_pfjets, true);
1065 TGraphErrors gr_pfjetsFwd = EresGraph(&plots_pfjets, false);
1066 gr_pfjets.SetName("pfJet");
1067 gr_pfjetsFwd.SetName("pfJetFwd");
1068
1069 // CaloJets Energy Resolution
1070 std::vector<resolPlot> plots_calojets;
1071 HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "CaloJet");
1072 GetJetsEres(&plots_calojets, branchCaloJet, branchGenJet, treeReaderJet);
1073 TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
1074 TGraphErrors gr_calojetsFwd = EresGraph(&plots_calojets, false);
1075 gr_calojets.SetName("caloJet");
1076 gr_calojetsFwd.SetName("caloJetFwd");
1077
1078 // MET Resolution vs HT
1079 std::vector<resolPlot> plots_met;
1080 HistogramsCollection(&plots_met, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "MET", -500, 500);
1081 GetMetres(&plots_met, branchScalarHT, branchMet, branchPFJet, treeReaderJet);
1082 TGraphErrors gr_met = EresGraph(&plots_met, true);
1083 gr_calojets.SetName("MET");
1084
1085 // Canvas
1086 TString btagEff = "btagEff";
1087 TCanvas *C_btag1 = new TCanvas(btagEff,btagEff, 1600, 600);
1088 C_btag1->Divide(2);
1089 C_btag1->cd(1);
1090 TLegend *leg_btag1 = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1091 leg_btag1->SetHeader("#splitline{B-tagging}{|#eta| < 1.5}");
1092 leg_btag1->AddEntry("","","");
1093
1094 gPad->SetLogx();
1095 histos_btag.first->Draw();
1096 addHist(histos_btag.first, leg_btag1, 0);
1097
1098 DrawAxis(histos_btag.first, leg_btag1, 0);
1099 leg_btag1->Draw();
1100
1101 C_btag1->cd(2);
1102 TLegend *leg_btag2 = new TLegend(resLegXmin,resLegYmin,resLegXmax+0.05,resLegYmax);
1103 leg_btag2->SetHeader("#splitline{B-tagging}{1.5 < |#eta| < 2.5}");
1104 leg_btag2->AddEntry("","","");
1105
1106 gPad->SetLogx();
1107 histos_btag.second->Draw();
1108 addHist(histos_btag.second, leg_btag2, 0);
1109
1110 DrawAxis(histos_btag.second, leg_btag2, 0);
1111 leg_btag2->Draw();
1112 C_btag1->cd(0);
1113
1114 TString tautagEff = "tautagEff";
1115 TCanvas *C_tautag1 = new TCanvas(tautagEff,tautagEff, 1600, 600);
1116 C_tautag1->Divide(2);
1117 C_tautag1->cd(1);
1118 TLegend *leg_tautag1 = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1119 leg_tautag1->SetHeader("#splitline{#tau-tagging}{|#eta| < 1.5}");
1120 leg_tautag1->AddEntry("","","");
1121
1122 gPad->SetLogx();
1123 histos_tautag.first->Draw();
1124 addHist(histos_tautag.first, leg_tautag1, 0);
1125
1126 DrawAxis(histos_tautag.first, leg_tautag1, 0);
1127 leg_tautag1->Draw();
1128
1129 C_tautag1->cd(2);
1130 TLegend *leg_tautag2 = new TLegend(resLegXmin,resLegYmin,resLegXmax+0.05,resLegYmax);
1131 leg_tautag2->SetHeader("#splitline{#tau-tagging}{1.5 < |#eta| < 2.5}");
1132 leg_tautag2->AddEntry("","","");
1133
1134 gPad->SetLogx();
1135 histos_tautag.second->Draw();
1136 addHist(histos_tautag.second, leg_tautag2, 0);
1137
1138 DrawAxis(histos_tautag.second, leg_tautag2, 0);
1139 leg_tautag2->Draw();
1140 C_tautag1->cd(0);
1141
1142 TString jetRes = "jetERes";
1143 TString jetResFwd = "jetEResFwd";
1144 TCanvas *C_jet = new TCanvas(jetRes,jetRes, 1600, 600);
1145 C_jet->Divide(2);
1146
1147 C_jet->cd(1);
1148 TMultiGraph *mg_jet = new TMultiGraph(jetRes,jetRes);
1149 TLegend *leg_jet = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1150 leg_jet->SetHeader("#splitline{jets}{|#eta| < 1.5}");
1151 leg_jet->AddEntry("","","");
1152
1153 addGraph(mg_jet, &gr_calojets, leg_jet, 3);
1154 addGraph(mg_jet, &gr_pfjets, leg_jet, 1);
1155
1156 mg_jet->Draw("ACX");
1157 leg_jet->Draw();
1158
1159 DrawAxis(mg_jet, leg_jet, 0.25);
1160
1161 C_jet->cd(2);
1162 TMultiGraph *mg_jetFwd = new TMultiGraph(jetResFwd,jetResFwd);
1163 TLegend *leg_jetFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1164 leg_jetFwd->SetHeader("#splitline{jets}{|#eta| < 1.5}");
1165 leg_jetFwd->AddEntry("","","");
1166
1167 addGraph(mg_jetFwd, &gr_calojetsFwd, leg_jetFwd, 3);
1168 addGraph(mg_jetFwd, &gr_pfjetsFwd, leg_jetFwd, 1);
1169
1170 mg_jetFwd->Draw("ACX");
1171 leg_jetFwd->Draw();
1172
1173 DrawAxis(mg_jetFwd, leg_jetFwd, 0.25);
1174
1175 TString metRes = "MetRes";
1176 TCanvas *C_met = new TCanvas(metRes,metRes, 800, 600);
1177
1178 TMultiGraph *mg_met = new TMultiGraph(metRes,metRes);
1179 TLegend *leg_met = new TLegend(topLeftLegXmin,topLeftLegYmin+0.2,topLeftLegXmax-0.2,topLeftLegYmax);
1180 leg_met->SetHeader("E_{T}^{miss}");
1181 leg_met->AddEntry("","","");
1182
1183
1184 addGraph(mg_met, &gr_met, leg_met, 0);
1185
1186 mg_met->Draw("ACX");
1187 leg_met->Draw();
1188
1189 DrawAxis(mg_met, leg_met, 300);
1190
1191 mg_met->GetXaxis()->SetTitle("H_{T} [GeV]");
1192 mg_met->GetYaxis()->SetTitle("#sigma(ME_{x}) [GeV]");
1193
1194 C_jet->Print("delphes_validation.pdf","pdf");
1195 C_btag1->Print("delphes_validation.pdf","pdf");
1196 C_tautag1->Print("delphes_validation.pdf","pdf");
1197 C_met->Print("delphes_validation.pdf)","pdf");
1198
1199 TFile *fout = new TFile(outputFile,"recreate");
1200
1201 for (int bin = 0; bin < Nbins; bin++)
1202 {
1203 plots_el.at(bin).cenResolHist->Write();
1204 plots_eltrack.at(bin).cenResolHist->Write();
1205 plots_eltower.at(bin).cenResolHist->Write();
1206 plots_el.at(bin).fwdResolHist->Write();
1207 plots_eltrack.at(bin).fwdResolHist->Write();
1208 plots_eltower.at(bin).fwdResolHist->Write();
1209 }
1210
1211 histos_el.first->Write();
1212 histos_el.second->Write();
1213 histos_eltrack.first->Write();
1214 histos_eltrack.second->Write();
1215 histos_eltower.first->Write();
1216 histos_eltower.second->Write();
1217 C_el1->Write();
1218 C_el2->Write();
1219
1220 histos_mu.first->Write();
1221 histos_mu.second->Write();
1222 histos_mutrack.first->Write();
1223 histos_mutrack.second->Write();
1224 C_mu1->Write();
1225 C_mu->Write();
1226
1227 histos_ph.first->Write();
1228 histos_ph.second->Write();
1229 C_ph1->Write();
1230 C_ph->Write();
1231
1232 for (int bin = 0; bin < Nbins; bin++)
1233 {
1234 plots_pfjets.at(bin).cenResolHist->Write();
1235 plots_pfjets.at(bin).fwdResolHist->Write();
1236 plots_calojets.at(bin).cenResolHist->Write();
1237 plots_calojets.at(bin).fwdResolHist->Write();
1238 plots_met.at(bin).cenResolHist->Write();
1239 }
1240 histos_btag.first->Write();
1241 histos_btag.second->Write();
1242 histos_tautag.first->Write();
1243 histos_tautag.second->Write();
1244 C_btag1->Write();
1245 C_tautag1->Write();
1246 C_jet->Write();
1247 C_met->Write();
1248
1249 fout->Write();
1250
1251 cout << "** Exiting..." << endl;
1252
1253 delete treeReaderElectron;
1254 delete treeReaderMuon;
1255 delete treeReaderPhoton;
1256 delete treeReaderJet;
1257 delete treeReaderBJet;
1258 delete treeReaderTauJet;
1259 delete chainElectron;
1260 delete chainMuon;
1261 delete chainPhoton;
1262 delete chainJet;
1263 delete chainBJet;
1264 delete chainTauJet;
1265}
1266
1267//------------------------------------------------------------------------------
1268
1269int main(int argc, char *argv[])
1270{
1271 char *appName = "Validation";
1272
1273 if(argc != 8)
1274 {
1275 cout << " Usage: " << appName << " input_file_electron input_file_muon input_file_photon input_file_jet input_file_bjet input_file_taujet output_file" << endl;
1276 cout << " input_file_electron - input file in ROOT format ('Delphes' tree)," << endl;
1277 cout << " input_file_muon - input file in ROOT format ('Delphes' tree)," << endl;
1278 cout << " input_file_photon - input file in ROOT format ('Delphes' tree)," << endl;
1279 cout << " input_file_jet - input file in ROOT format ('Delphes' tree)," << endl;
1280 cout << " input_file_bjet - input file in ROOT format ('Delphes' tree)," << endl;
1281 cout << " input_file_taujet - input file in ROOT format ('Delphes' tree)," << endl;
1282 cout << " output_file - output file in ROOT format" << endl;
1283 return 1;
1284 }
1285
1286 gROOT->SetBatch();
1287
1288 int appargc = 1;
1289 char *appargv[] = {appName};
1290 TApplication app(appName, &appargc, appargv);
1291
1292 Validation(argv[1], argv[2], argv[3], argv[4], argv[5], argv[6], argv[7]);
1293}
1294
1295
Note: See TracBrowser for help on using the repository browser.