Fork me on GitHub

source: git/examples/Validation.cpp@ a17ec5a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a17ec5a was 360d17e, checked in by Pavel Demin <pavel.demin@…>, 8 years ago

improve validation code

  • Property mode set to 100644
File size: 38.9 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 Double_t sig = f1->GetParameter(2);
479 Double_t sigErr = f1->GetParError(2);
480 delete f1;
481 return make_pair (sig, sigErr);
482}
483
484
485TGraphErrors EresGraph(std::vector<resolPlot> *histos, bool central, bool rms = false)
486{
487 Int_t bin;
488 Int_t count = 0;
489 TGraphErrors gr = TGraphErrors(Nbins/2);
490 Double_t sig = 0;
491 Double_t sigErr = 0;
492 for (bin = 0; bin < Nbins; bin++)
493 {
494 if (central == true && histos->at(bin).cenResolHist->GetEntries() > 100)
495 {
496 std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).cenResolHist);
497 if (rms == true)
498 {
499 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
500 gr.SetPointError(count,0, sigvalues.second); // to correct
501 }
502 else
503 {
504 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
505 gr.SetPointError(count,0, sigvalues.second);
506 }
507 count++;
508 }
509
510 else if (central == false && histos->at(bin).fwdResolHist->GetEntries() > 10)
511 {
512 std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).fwdResolHist);
513 if (rms == true)
514 {
515 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
516 gr.SetPointError(count,0, sigvalues.second); // to correct
517 }
518 else
519 {
520 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
521 gr.SetPointError(count,0, sigvalues.second);
522 }
523 count++;
524 }
525
526 }
527 return gr;
528}
529
530
531//------------------------------------------------------------------------------
532
533
534// type 1 : object, 2 : track, 3 : tower
535
536void addGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int type)
537{
538
539 gr->SetLineWidth(2);
540
541 switch ( type )
542 {
543 case 1:
544 gr->SetLineColor(objColor);
545 gr->SetLineStyle(objStyle);
546 std::cout << "Adding " << gr->GetName() << std::endl;
547 mg->Add(gr);
548 leg->AddEntry(gr,"Reco","l");
549 break;
550
551 case 2:
552 gr->SetLineColor(trackColor);
553 gr->SetLineStyle(trackStyle);
554 mg->Add(gr);
555 leg->AddEntry(gr,"Track","l");
556 break;
557
558 case 3:
559 gr->SetLineColor(towerColor);
560 gr->SetLineStyle(towerStyle);
561 mg->Add(gr);
562 leg->AddEntry(gr,"Tower","l");
563 break;
564
565 case 0:
566 gr->SetLineColor(objColor);
567 gr->SetLineStyle(objStyle);
568 mg->Add(gr);
569 break;
570
571 default:
572 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
573 break;
574 }
575}
576
577void addHist(TH1D *h, TLegend *leg, int type)
578{
579 h->SetLineWidth(2);
580
581 switch ( type )
582 {
583 case 1:
584 h->SetLineColor(objColor);
585 h->SetLineStyle(objStyle);
586 leg->AddEntry(h,"Reco","l");
587 break;
588
589 case 2:
590 h->SetLineColor(trackColor);
591 h->SetLineStyle(trackStyle);
592 leg->AddEntry(h,"Track","l");
593 break;
594
595 case 3:
596 h->SetLineColor(towerColor);
597 h->SetLineStyle(towerStyle);
598 leg->AddEntry(h,"Tower","l");
599 break;
600
601 case 0:
602 h->SetLineColor(objColor);
603 h->SetLineStyle(objStyle);
604 break;
605
606 default:
607 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
608 break;
609 }
610}
611
612void DrawAxis(TMultiGraph *mg, TLegend *leg, double max)
613{
614 mg->SetMinimum(0.);
615 mg->SetMaximum(max);
616 mg->GetXaxis()->SetLimits(ptrangemin,ptrangemax);
617 mg->GetYaxis()->SetTitle("resolution");
618 mg->GetXaxis()->SetTitle("p_{T} [GeV]");
619 mg->GetYaxis()->SetTitleSize(0.07);
620 mg->GetXaxis()->SetTitleSize(0.07);
621 mg->GetYaxis()->SetLabelSize(0.06);
622 mg->GetXaxis()->SetLabelSize(0.06);
623 mg->GetYaxis()->SetLabelOffset(0.03);
624 mg->GetYaxis()->SetTitleOffset(1.4);
625 mg->GetXaxis()->SetTitleOffset(1.4);
626
627 mg->GetYaxis()->SetNdivisions(505);
628
629 leg->SetBorderSize(0);
630 leg->SetShadowColor(0);
631 leg->SetFillColor(0);
632 leg->SetFillStyle(0);
633
634 gStyle->SetOptTitle(0);
635 gPad->SetLogx();
636 gPad->SetBottomMargin(0.2);
637 gPad->SetLeftMargin(0.2);
638 gPad->Modified();
639 gPad->Update();
640
641}
642
643void DrawAxis(TH1D *h, TLegend *leg, int type)
644{
645
646 h->GetYaxis()->SetRangeUser(0,1.0);
647 if (type == 0) h->GetXaxis()->SetTitle("p_{T} [GeV]");
648 else h->GetXaxis()->SetTitle("#eta");
649 h->GetYaxis()->SetTitle("efficiency");
650 h->GetYaxis()->SetTitleSize(0.07);
651 h->GetXaxis()->SetTitleSize(0.07);
652 h->GetYaxis()->SetLabelSize(0.06);
653 h->GetXaxis()->SetLabelSize(0.06);
654 h->GetYaxis()->SetLabelOffset(0.03);
655 h->GetYaxis()->SetTitleOffset(1.3);
656 h->GetXaxis()->SetTitleOffset(1.4);
657
658 h->GetYaxis()->SetNdivisions(505);
659
660 leg->SetBorderSize(0);
661 leg->SetShadowColor(0);
662 leg->SetFillColor(0);
663 leg->SetFillStyle(0);
664
665 gStyle->SetOptTitle(0);
666 gStyle->SetOptStat(0);
667 gPad->SetBottomMargin(0.2);
668 gPad->SetLeftMargin(0.2);
669
670 gPad->Modified();
671 gPad->Update();
672
673}
674
675
676void Validation(const char *inputFileElectron, const char *inputFileMuon, const char *inputFilePhoton, const char *inputFileJet, const char *inputFileBJet, const char *inputFileTauJet, const char *outputFile)
677{
678 TChain *chainElectron = new TChain("Delphes");
679 chainElectron->Add(inputFileElectron);
680 ExRootTreeReader *treeReaderElectron = new ExRootTreeReader(chainElectron);
681
682 TChain *chainMuon = new TChain("Delphes");
683 chainMuon->Add(inputFileMuon);
684 ExRootTreeReader *treeReaderMuon = new ExRootTreeReader(chainMuon);
685
686 TChain *chainPhoton = new TChain("Delphes");
687 chainPhoton->Add(inputFilePhoton);
688 ExRootTreeReader *treeReaderPhoton = new ExRootTreeReader(chainPhoton);
689
690 TChain *chainJet = new TChain("Delphes");
691 chainJet->Add(inputFileJet);
692 ExRootTreeReader *treeReaderJet = new ExRootTreeReader(chainJet);
693
694 TChain *chainBJet = new TChain("Delphes");
695 chainBJet->Add(inputFileBJet);
696 ExRootTreeReader *treeReaderBJet = new ExRootTreeReader(chainBJet);
697
698 TChain *chainTauJet = new TChain("Delphes");
699 chainTauJet->Add(inputFileTauJet);
700 ExRootTreeReader *treeReaderTauJet = new ExRootTreeReader(chainTauJet);
701
702 TClonesArray *branchParticleElectron = treeReaderElectron->UseBranch("Particle");
703 TClonesArray *branchTrackElectron = treeReaderElectron->UseBranch("Track");
704 TClonesArray *branchTowerElectron = treeReaderElectron->UseBranch("Tower");
705 TClonesArray *branchElectron = treeReaderElectron->UseBranch("Electron");
706
707 TClonesArray *branchParticleMuon = treeReaderMuon->UseBranch("Particle");
708 TClonesArray *branchTrackMuon = treeReaderMuon->UseBranch("Track");
709 TClonesArray *branchMuon = treeReaderMuon->UseBranch("Muon");
710
711 TClonesArray *branchParticlePhoton = treeReaderPhoton->UseBranch("Particle");
712 TClonesArray *branchTowerPhoton = treeReaderPhoton->UseBranch("Tower");
713 TClonesArray *branchPhoton = treeReaderPhoton->UseBranch("Photon");
714
715 TClonesArray *branchGenJet = treeReaderJet->UseBranch("GenJet");
716 TClonesArray *branchPFJet = treeReaderJet->UseBranch("Jet");
717 TClonesArray *branchCaloJet = treeReaderJet->UseBranch("CaloJet");
718
719 TClonesArray *branchParticleBJet = treeReaderBJet->UseBranch("Particle");
720 TClonesArray *branchPFBJet = treeReaderBJet->UseBranch("Jet");
721
722 TClonesArray *branchParticleTauJet = treeReaderTauJet->UseBranch("Particle");
723 TClonesArray *branchPFTauJet = treeReaderTauJet->UseBranch("Jet");
724
725 TClonesArray *branchScalarHT = treeReaderJet->UseBranch("ScalarHT");
726 TClonesArray *branchMet = treeReaderJet->UseBranch("MissingET");
727
728 ///////////////
729 // Electrons //
730 ///////////////
731
732 // Reconstruction efficiency
733 TString elecs = "Electron";
734 int elID = 11;
735 std::pair<TH1D*,TH1D*> histos_el = GetEff<Electron>(branchElectron, branchParticleElectron, "Electron", elID, treeReaderElectron);
736
737 // tracking reconstruction efficiency
738 std::pair <TH1D*,TH1D*> histos_eltrack = GetEff<Track>(branchTrackElectron, branchParticleElectron, "electronTrack", elID, treeReaderElectron);
739
740 // Tower reconstruction efficiency
741 std::pair <TH1D*,TH1D*> histos_eltower = GetEff<Tower>(branchTowerElectron, branchParticleElectron, "electronTower", elID, treeReaderElectron);
742
743 // Electron Energy Resolution
744 std::vector<resolPlot> plots_el;
745 HistogramsCollection(&plots_el, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electrons");
746 GetEres<Electron>(&plots_el, branchElectron, branchParticleElectron, elID, treeReaderElectron);
747 TGraphErrors gr_el = EresGraph(&plots_el, true);
748 TGraphErrors gr_elFwd = EresGraph(&plots_el, false);
749 gr_el.SetName("Electron");
750 gr_elFwd.SetName("ElectronFwd");
751
752 // Electron Track Energy Resolution
753 std::vector<resolPlot> plots_eltrack;
754 HistogramsCollection(&plots_eltrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTracks");
755 GetEres<Track>(&plots_eltrack, branchTrackElectron, branchParticleElectron, elID, treeReaderElectron);
756 TGraphErrors gr_eltrack = EresGraph(&plots_eltrack, true);
757 TGraphErrors gr_eltrackFwd = EresGraph(&plots_eltrack, false);
758 gr_eltrack.SetName("ElectronTracks");
759 gr_eltrackFwd.SetName("ElectronTracksFwd");
760
761 // Electron Tower Energy Resolution
762 std::vector<resolPlot> plots_eltower;
763 HistogramsCollection(&plots_eltower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTower");
764 GetEres<Tower>(&plots_eltower, branchTowerElectron, branchParticleElectron, elID, treeReaderElectron);
765 TGraphErrors gr_eltower = EresGraph(&plots_eltower, true);
766 TGraphErrors gr_eltowerFwd = EresGraph(&plots_eltower, false);
767 gr_eltower.SetName("ElectronTower");
768 gr_eltrackFwd.SetName("ElectronTracksFwd");
769
770 // Canvases
771 TString elEff = "electronEff";
772 TCanvas *C_el1 = new TCanvas(elEff,elEff, 1600, 600);
773 C_el1->Divide(2);
774 C_el1->cd(1);
775 TLegend *leg_el1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
776 leg_el1->SetHeader("#splitline{electrons}{|#eta| < 1.5}");
777 leg_el1->AddEntry("","","");
778
779 gPad->SetLogx();
780 histos_eltrack.first->Draw("][");
781 addHist(histos_eltrack.first, leg_el1, 2);
782 histos_el.first->Draw("same ][");
783 addHist(histos_el.first, leg_el1, 1);
784 DrawAxis(histos_eltrack.first, leg_el1, 0);
785
786 leg_el1->Draw();
787
788 C_el1->cd(2);
789 TLegend *leg_el2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
790 leg_el2->SetHeader("#splitline{electrons}{1.5 < |#eta| < 2.5}");
791 leg_el2->AddEntry("","","");
792
793 gPad->SetLogx();
794 histos_eltrack.second->Draw("][");
795 addHist(histos_eltrack.second, leg_el2, 2);
796 histos_el.second->Draw("same ][");
797 addHist(histos_el.second, leg_el2, 1);
798
799 DrawAxis(histos_eltrack.second, leg_el2, 0);
800 leg_el2->Draw();
801
802 TString elRes = "electronERes";
803 TString elResFwd = "electronEResForward";
804 TCanvas *C_el2 = new TCanvas(elRes,elRes, 1600, 600);
805 C_el2->Divide(2);
806 C_el2->cd(1);
807 TMultiGraph *mg_el = new TMultiGraph(elRes,elRes);
808 TLegend *leg_el = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
809 leg_el->SetHeader("#splitline{electrons}{|#eta| < 1.5}");
810 leg_el->AddEntry("","","");
811
812 addGraph(mg_el, &gr_eltower, leg_el, 3);
813 addGraph(mg_el, &gr_eltrack, leg_el, 2);
814 addGraph(mg_el, &gr_el, leg_el, 1);
815
816 mg_el->Draw("ACX");
817 leg_el->Draw();
818
819 DrawAxis(mg_el, leg_el, 0.1);
820
821 C_el2->cd(2);
822 TMultiGraph *mg_elFwd = new TMultiGraph(elResFwd,elResFwd);
823 TLegend *leg_elFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
824 leg_elFwd->SetHeader("#splitline{electrons}{1.5 < |#eta| < 2.5}");
825 leg_elFwd->AddEntry("","","");
826
827 addGraph(mg_elFwd, &gr_eltowerFwd, leg_elFwd, 3);
828 addGraph(mg_elFwd, &gr_eltrackFwd, leg_elFwd, 2);
829 addGraph(mg_elFwd, &gr_elFwd, leg_elFwd, 1);
830
831 mg_elFwd->Draw("ACX");
832 leg_elFwd->Draw();
833
834 DrawAxis(mg_elFwd, leg_elFwd, 0.2);
835
836 C_el1->Print("delphes_validation.pdf(","pdf");
837 C_el2->Print("delphes_validation.pdf","pdf");
838
839 gDirectory->cd(0);
840
841 ///////////
842 // Muons //
843 ///////////
844
845 // Reconstruction efficiency
846 int muID = 13;
847 std::pair<TH1D*,TH1D*> histos_mu = GetEff<Muon>(branchMuon, branchParticleMuon,"Muon", muID, treeReaderMuon);
848
849 // muon tracking reconstruction efficiency
850 std::pair <TH1D*,TH1D*> histos_mutrack = GetEff<Track>(branchTrackMuon, branchParticleMuon, "muonTrack", muID, treeReaderMuon);
851
852 // Muon Energy Resolution
853 std::vector<resolPlot> plots_mu;
854 HistogramsCollection(&plots_mu, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muons");
855 GetEres<Muon>(&plots_mu, branchMuon, branchParticleMuon, muID, treeReaderMuon);
856 TGraphErrors gr_mu = EresGraph(&plots_mu, true);
857 TGraphErrors gr_muFwd = EresGraph(&plots_mu, false);
858 gr_mu.SetName("Muon");
859 gr_muFwd.SetName("MuonFwd");
860
861 // Muon Track Energy Resolution
862 std::vector<resolPlot> plots_mutrack;
863 HistogramsCollection(&plots_mutrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muonsTracks");
864 GetEres<Track>(&plots_mutrack, branchTrackMuon, branchParticleMuon, muID, treeReaderMuon);
865 TGraphErrors gr_mutrack = EresGraph(&plots_mutrack, true);
866 TGraphErrors gr_mutrackFwd = EresGraph(&plots_mutrack, false);
867 gr_mutrackFwd.SetName("MuonTracksFwd");
868
869 // Canvas
870
871 TString muEff = "muonEff";
872 TCanvas *C_mu1 = new TCanvas(muEff,muEff, 1600, 600);
873 C_mu1->Divide(2);
874 C_mu1->cd(1);
875 TLegend *leg_mu1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
876 leg_mu1->SetHeader("#splitline{muons}{|#eta| < 1.5}");
877 leg_mu1->AddEntry("","","");
878
879
880 gPad->SetLogx();
881 histos_mutrack.first->Draw("][");
882 addHist(histos_mutrack.first, leg_mu1, 2);
883 histos_mu.first->Draw("same ][");
884 addHist(histos_mu.first, leg_mu1, 1);
885
886 DrawAxis(histos_mutrack.first, leg_mu1, 0);
887
888 leg_mu1->Draw();
889
890 C_mu1->cd(2);
891 TLegend *leg_mu2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
892 leg_mu2->SetHeader("#splitline{muons}{1.5 < |#eta| < 2.5}");
893 leg_mu2->AddEntry("","","");
894
895 gPad->SetLogx();
896 histos_mutrack.second->Draw("][");
897 addHist(histos_mutrack.second, leg_mu2, 2);
898 histos_mu.second->Draw("same ][");
899 addHist(histos_mu.second, leg_mu2, 1);
900
901 DrawAxis(histos_mutrack.second, leg_mu2, 1);
902 leg_mu2->Draw();
903
904 TString muRes = "muonERes";
905 TString muResFwd = "muonEResFwd";
906
907 TCanvas *C_mu = new TCanvas(muRes,muRes, 1600, 600);
908 C_mu->Divide(2);
909 C_mu->cd(1);
910 TMultiGraph *mg_mu = new TMultiGraph(muRes,muRes);
911 TLegend *leg_mu = new TLegend(topLeftLegXmin,topLeftLegYmin,topLeftLegXmax,topLeftLegYmax);
912 leg_mu->SetHeader("#splitline{muons}{|#eta| < 1.5}");
913 leg_mu->AddEntry("","","");
914
915 addGraph(mg_mu, &gr_mutrack, leg_mu, 2);
916 addGraph(mg_mu, &gr_mu, leg_mu, 1);
917
918 mg_mu->Draw("ACX");
919 leg_mu->Draw();
920
921 DrawAxis(mg_mu, leg_mu, 0.3);
922
923 C_mu->cd(2);
924 TMultiGraph *mg_muFwd = new TMultiGraph(muResFwd,muResFwd);
925 TLegend *leg_muFwd = new TLegend(topLeftLegXmin,topLeftLegYmin,topLeftLegXmax,topLeftLegYmax);
926 leg_muFwd->SetHeader("#splitline{muons}{1.5 < |#eta| < 2.5}");
927 leg_muFwd->AddEntry("","","");
928
929 addGraph(mg_muFwd, &gr_mutrackFwd, leg_muFwd, 2);
930 addGraph(mg_muFwd, &gr_muFwd, leg_muFwd, 1);
931
932 mg_muFwd->Draw("ACX");
933 leg_muFwd->Draw();
934
935 DrawAxis(mg_muFwd, leg_muFwd, 0.3);
936
937 C_mu1->Print("delphes_validation.pdf","pdf");
938 C_mu->Print("delphes_validation.pdf","pdf");
939
940 gDirectory->cd(0);
941
942 /////////////
943 // Photons //
944 /////////////
945
946 // Reconstruction efficiency
947 int phID = 22;
948 std::pair<TH1D*,TH1D*> histos_ph = GetEff<Electron>(branchPhoton, branchParticlePhoton, "Photon", phID, treeReaderPhoton);
949 std::pair<TH1D*,TH1D*> histos_phtower = GetEff<Electron>(branchTowerPhoton, branchParticlePhoton, "Photon", phID, treeReaderPhoton);
950
951 // Photon Energy Resolution
952 std::vector<resolPlot> plots_ph;
953 HistogramsCollection(&plots_ph, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photons");
954 GetEres<Photon>(&plots_ph, branchPhoton, branchParticlePhoton, phID, treeReaderPhoton);
955 TGraphErrors gr_ph = EresGraph(&plots_ph, true);
956 TGraphErrors gr_phFwd = EresGraph(&plots_ph, false);
957 gr_ph.SetName("Photon");
958 gr_phFwd.SetName("PhotonFwd");
959
960
961 // Photon Tower Energy Resolution
962 std::vector<resolPlot> plots_phtower;
963 HistogramsCollection(&plots_phtower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photonsTower");
964 GetEres<Tower>(&plots_phtower, branchTowerPhoton, branchParticlePhoton, phID, treeReaderPhoton);
965 TGraphErrors gr_phtower = EresGraph(&plots_phtower, true);
966 TGraphErrors gr_phtowerFwd = EresGraph(&plots_phtower, false);
967 gr_phtower.SetName("PhotonTower");
968 gr_phtowerFwd.SetName("PhotonTowerFwd");
969
970 // Canvas
971
972 TString phEff = "photonEff";
973 TCanvas *C_ph1 = new TCanvas(phEff,phEff, 1600, 600);
974 C_ph1->Divide(2);
975 C_ph1->cd(1);
976 TLegend *leg_ph1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
977 leg_ph1->SetHeader("#splitline{photons}{|#eta| < 1.5}");
978 leg_ph1->AddEntry("","","");
979
980
981 gPad->SetLogx();
982 histos_phtower.first->Draw("][");
983 addHist(histos_phtower.first, leg_ph1, 3);
984 histos_ph.first->Draw("same ][");
985 addHist(histos_ph.first, leg_ph1, 1);
986
987 DrawAxis(histos_phtower.first, leg_ph1, 0);
988 leg_ph1->Draw();
989
990 C_ph1->cd(2);
991 TLegend *leg_ph2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
992 leg_ph2->SetHeader("#splitline{photons}{1.5 < |#eta| < 2.5}");
993 leg_ph2->AddEntry("","","");
994
995
996 gPad->SetLogx();
997 histos_phtower.second->Draw("][");
998 addHist(histos_phtower.second, leg_ph2, 3);
999 histos_ph.second->Draw("same ][");
1000 addHist(histos_ph.second, leg_ph2, 1);
1001
1002 DrawAxis(histos_phtower.second, leg_ph2, 1);
1003 leg_ph2->Draw();
1004
1005 TString phRes = "phERes";
1006 TString phResFwd = "phEResFwd";
1007
1008 TCanvas *C_ph = new TCanvas(phRes,phRes, 1600, 600);
1009 C_ph->Divide(2);
1010 C_ph->cd(1);
1011 TMultiGraph *mg_ph = new TMultiGraph(phRes,phRes);
1012 TLegend *leg_ph = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1013 leg_ph->SetHeader("#splitline{photons}{|#eta| < 1.5}");
1014 leg_ph->AddEntry("","","");
1015
1016 addGraph(mg_ph, &gr_phtower, leg_ph, 3);
1017 addGraph(mg_ph, &gr_ph, leg_ph, 1);
1018
1019 mg_ph->Draw("ACX");
1020 leg_ph->Draw();
1021
1022 DrawAxis(mg_ph, leg_ph, 0.3);
1023
1024 C_ph->cd(2);
1025 TMultiGraph *mg_phFwd = new TMultiGraph(phResFwd,phResFwd);
1026 TLegend *leg_phFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1027 leg_phFwd->SetHeader("#splitline{photons}{1.5 < |#eta| < 2.5}");
1028 leg_phFwd->AddEntry("","","");
1029
1030 addGraph(mg_phFwd, &gr_phtowerFwd, leg_phFwd, 3);
1031 addGraph(mg_phFwd, &gr_phFwd, leg_phFwd, 1);
1032
1033 mg_phFwd->Draw("ACX");
1034 leg_phFwd->Draw();
1035
1036 DrawAxis(mg_phFwd, leg_phFwd, 0.3);
1037
1038 C_ph1->Print("delphes_validation.pdf","pdf");
1039 C_ph->Print("delphes_validation.pdf","pdf");
1040
1041 gDirectory->cd(0);
1042
1043 //////////
1044 // Jets //
1045 //////////
1046
1047 // BJets Reconstruction efficiency
1048 int bID = 5;
1049 std::pair<TH1D*,TH1D*> histos_btag = GetEff<Jet>(branchPFBJet, branchParticleBJet,"BTag", bID, treeReaderBJet);
1050
1051 // TauJets Reconstruction efficiency
1052 int tauID = 15;
1053 std::pair<TH1D*,TH1D*> histos_tautag = GetEff<Jet>(branchPFTauJet, branchParticleTauJet,"TauTag", tauID, treeReaderTauJet);
1054
1055 // PFJets Energy Resolution
1056 std::vector<resolPlot> plots_pfjets;
1057 HistogramsCollection(&plots_pfjets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "PFJet");
1058 GetJetsEres(&plots_pfjets, branchPFJet, branchGenJet, treeReaderJet);
1059 TGraphErrors gr_pfjets = EresGraph(&plots_pfjets, true);
1060 TGraphErrors gr_pfjetsFwd = EresGraph(&plots_pfjets, false);
1061 gr_pfjets.SetName("pfJet");
1062 gr_pfjetsFwd.SetName("pfJetFwd");
1063
1064 // CaloJets Energy Resolution
1065 std::vector<resolPlot> plots_calojets;
1066 HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "CaloJet");
1067 GetJetsEres(&plots_calojets, branchCaloJet, branchGenJet, treeReaderJet);
1068 TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
1069 TGraphErrors gr_calojetsFwd = EresGraph(&plots_calojets, false);
1070 gr_calojets.SetName("caloJet");
1071 gr_calojetsFwd.SetName("caloJetFwd");
1072
1073 // MET Resolution vs HT
1074 std::vector<resolPlot> plots_met;
1075 HistogramsCollection(&plots_met, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "MET", -500, 500);
1076 GetMetres(&plots_met, branchScalarHT, branchMet, branchPFJet, treeReaderJet);
1077 TGraphErrors gr_met = EresGraph(&plots_met, true);
1078 gr_calojets.SetName("MET");
1079
1080 // Canvas
1081 TString btagEff = "btagEff";
1082 TCanvas *C_btag1 = new TCanvas(btagEff,btagEff, 1600, 600);
1083 C_btag1->Divide(2);
1084 C_btag1->cd(1);
1085 TLegend *leg_btag1 = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1086 leg_btag1->SetHeader("#splitline{B-tagging}{|#eta| < 1.5}");
1087 leg_btag1->AddEntry("","","");
1088
1089 gPad->SetLogx();
1090 histos_btag.first->Draw();
1091 addHist(histos_btag.first, leg_btag1, 0);
1092
1093 DrawAxis(histos_btag.first, leg_btag1, 0);
1094 leg_btag1->Draw();
1095
1096 C_btag1->cd(2);
1097 TLegend *leg_btag2 = new TLegend(resLegXmin,resLegYmin,resLegXmax+0.05,resLegYmax);
1098 leg_btag2->SetHeader("#splitline{B-tagging}{1.5 < |#eta| < 2.5}");
1099 leg_btag2->AddEntry("","","");
1100
1101 gPad->SetLogx();
1102 histos_btag.second->Draw();
1103 addHist(histos_btag.second, leg_btag2, 0);
1104
1105 DrawAxis(histos_btag.second, leg_btag2, 0);
1106 leg_btag2->Draw();
1107 C_btag1->cd(0);
1108
1109 TString tautagEff = "tautagEff";
1110 TCanvas *C_tautag1 = new TCanvas(tautagEff,tautagEff, 1600, 600);
1111 C_tautag1->Divide(2);
1112 C_tautag1->cd(1);
1113 TLegend *leg_tautag1 = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1114 leg_tautag1->SetHeader("#splitline{#tau-tagging}{|#eta| < 1.5}");
1115 leg_tautag1->AddEntry("","","");
1116
1117 gPad->SetLogx();
1118 histos_tautag.first->Draw();
1119 addHist(histos_tautag.first, leg_tautag1, 0);
1120
1121 DrawAxis(histos_tautag.first, leg_tautag1, 0);
1122 leg_tautag1->Draw();
1123
1124 C_tautag1->cd(2);
1125 TLegend *leg_tautag2 = new TLegend(resLegXmin,resLegYmin,resLegXmax+0.05,resLegYmax);
1126 leg_tautag2->SetHeader("#splitline{#tau-tagging}{1.5 < |#eta| < 2.5}");
1127 leg_tautag2->AddEntry("","","");
1128
1129 gPad->SetLogx();
1130 histos_tautag.second->Draw();
1131 addHist(histos_tautag.second, leg_tautag2, 0);
1132
1133 DrawAxis(histos_tautag.second, leg_tautag2, 0);
1134 leg_tautag2->Draw();
1135 C_tautag1->cd(0);
1136
1137 TString jetRes = "jetERes";
1138 TString jetResFwd = "jetEResFwd";
1139 TCanvas *C_jet = new TCanvas(jetRes,jetRes, 1600, 600);
1140 C_jet->Divide(2);
1141
1142 C_jet->cd(1);
1143 TMultiGraph *mg_jet = new TMultiGraph(jetRes,jetRes);
1144 TLegend *leg_jet = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1145 leg_jet->SetHeader("#splitline{jets}{|#eta| < 1.5}");
1146 leg_jet->AddEntry("","","");
1147
1148 addGraph(mg_jet, &gr_calojets, leg_jet, 3);
1149 addGraph(mg_jet, &gr_pfjets, leg_jet, 1);
1150
1151 mg_jet->Draw("ACX");
1152 leg_jet->Draw();
1153
1154 DrawAxis(mg_jet, leg_jet, 0.25);
1155
1156 C_jet->cd(2);
1157 TMultiGraph *mg_jetFwd = new TMultiGraph(jetResFwd,jetResFwd);
1158 TLegend *leg_jetFwd = new TLegend(resLegXmin,resLegYmin,resLegXmax,resLegYmax);
1159 leg_jetFwd->SetHeader("#splitline{jets}{|#eta| < 1.5}");
1160 leg_jetFwd->AddEntry("","","");
1161
1162 addGraph(mg_jetFwd, &gr_calojetsFwd, leg_jetFwd, 3);
1163 addGraph(mg_jetFwd, &gr_pfjetsFwd, leg_jetFwd, 1);
1164
1165 mg_jetFwd->Draw("ACX");
1166 leg_jetFwd->Draw();
1167
1168 DrawAxis(mg_jetFwd, leg_jetFwd, 0.25);
1169
1170 TString metRes = "MetRes";
1171 TCanvas *C_met = new TCanvas(metRes,metRes, 800, 600);
1172
1173 TMultiGraph *mg_met = new TMultiGraph(metRes,metRes);
1174 TLegend *leg_met = new TLegend(topLeftLegXmin,topLeftLegYmin+0.2,topLeftLegXmax-0.2,topLeftLegYmax);
1175 leg_met->SetHeader("E_{T}^{miss}");
1176 leg_met->AddEntry("","","");
1177
1178
1179 addGraph(mg_met, &gr_met, leg_met, 0);
1180
1181 mg_met->Draw("ACX");
1182 leg_met->Draw();
1183
1184 DrawAxis(mg_met, leg_met, 300);
1185
1186 mg_met->GetXaxis()->SetTitle("H_{T} [GeV]");
1187 mg_met->GetYaxis()->SetTitle("#sigma(ME_{x}) [GeV]");
1188
1189 C_jet->Print("delphes_validation.pdf","pdf");
1190 C_btag1->Print("delphes_validation.pdf","pdf");
1191 C_tautag1->Print("delphes_validation.pdf","pdf");
1192 C_met->Print("delphes_validation.pdf)","pdf");
1193
1194 TFile *fout = new TFile(outputFile,"recreate");
1195
1196 for (int bin = 0; bin < Nbins; bin++)
1197 {
1198 plots_el.at(bin).cenResolHist->Write();
1199 plots_eltrack.at(bin).cenResolHist->Write();
1200 plots_eltower.at(bin).cenResolHist->Write();
1201 plots_el.at(bin).fwdResolHist->Write();
1202 plots_eltrack.at(bin).fwdResolHist->Write();
1203 plots_eltower.at(bin).fwdResolHist->Write();
1204 }
1205
1206 histos_el.first->Write();
1207 histos_el.second->Write();
1208 histos_eltrack.first->Write();
1209 histos_eltrack.second->Write();
1210 histos_eltower.first->Write();
1211 histos_eltower.second->Write();
1212 C_el1->Write();
1213 C_el2->Write();
1214
1215 histos_mu.first->Write();
1216 histos_mu.second->Write();
1217 histos_mutrack.first->Write();
1218 histos_mutrack.second->Write();
1219 C_mu1->Write();
1220 C_mu->Write();
1221
1222 histos_ph.first->Write();
1223 histos_ph.second->Write();
1224 C_ph1->Write();
1225 C_ph->Write();
1226
1227 for (int bin = 0; bin < Nbins; bin++)
1228 {
1229 plots_pfjets.at(bin).cenResolHist->Write();
1230 plots_pfjets.at(bin).fwdResolHist->Write();
1231 plots_calojets.at(bin).cenResolHist->Write();
1232 plots_calojets.at(bin).fwdResolHist->Write();
1233 plots_met.at(bin).cenResolHist->Write();
1234 }
1235 histos_btag.first->Write();
1236 histos_btag.second->Write();
1237 histos_tautag.first->Write();
1238 histos_tautag.second->Write();
1239 C_btag1->Write();
1240 C_tautag1->Write();
1241 C_jet->Write();
1242 C_met->Write();
1243
1244 fout->Write();
1245
1246 cout << "** Exiting..." << endl;
1247
1248 delete treeReaderElectron;
1249 delete treeReaderMuon;
1250 delete treeReaderPhoton;
1251 delete treeReaderJet;
1252 delete treeReaderBJet;
1253 delete treeReaderTauJet;
1254 delete chainElectron;
1255 delete chainMuon;
1256 delete chainPhoton;
1257 delete chainJet;
1258 delete chainBJet;
1259 delete chainTauJet;
1260}
1261
1262//------------------------------------------------------------------------------
1263
1264int main(int argc, char *argv[])
1265{
1266 char *appName = "Validation";
1267
1268 if(argc != 8)
1269 {
1270 cout << " Usage: " << appName << " input_file_electron input_file_muon input_file_photon input_file_jet input_file_bjet input_file_taujet output_file" << endl;
1271 cout << " input_file_electron - input file in ROOT format ('Delphes' tree)," << endl;
1272 cout << " input_file_muon - input file in ROOT format ('Delphes' tree)," << endl;
1273 cout << " input_file_photon - input file in ROOT format ('Delphes' tree)," << endl;
1274 cout << " input_file_jet - input file in ROOT format ('Delphes' tree)," << endl;
1275 cout << " input_file_bjet - input file in ROOT format ('Delphes' tree)," << endl;
1276 cout << " input_file_taujet - input file in ROOT format ('Delphes' tree)," << endl;
1277 cout << " output_file - output file in ROOT format" << endl;
1278 return 1;
1279 }
1280
1281 gROOT->SetBatch();
1282
1283 int appargc = 1;
1284 char *appargv[] = {appName};
1285 TApplication app(appName, &appargc, appargv);
1286
1287 Validation(argv[1], argv[2], argv[3], argv[4], argv[5], argv[6], argv[7]);
1288}
1289
1290
Note: See TracBrowser for help on using the repository browser.