Fork me on GitHub

source: git/examples/Validation.C@ 6563d4d

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

added met plot

  • Property mode set to 100644
File size: 29.6 KB
Line 
1#ifdef __CLING__
2R__LOAD_LIBRARY(libDelphes)
3#include "classes/DelphesClasses.h"
4#include "external/ExRootAnalysis/ExRootTreeReader.h"
5#include "external/ExRootAnalysis/ExRootResult.h"
6#else
7class ExRootTreeReader;
8class ExRootResult;
9#endif
10
11#include "TCanvas.h"
12#include "TSystem.h"
13#include "TStyle.h"
14#include "TLegend.h"
15#include <TH1.h>
16#include "TString.h"
17#include "vector"
18#include <TMath.h>
19#include <iostream>
20#include "TGraph.h"
21#include "TGraphErrors.h"
22#include "TMultiGraph.h"
23#include <typeinfo>
24#include "TLorentzVector.h"
25
26//------------------------------------------------------------------------------
27
28double ptrangemin = 10;
29double ptrangemax = 10000;
30static const int Nbins = 20;
31
32int objStyle = 1;
33int trackStyle = 7;
34int towerStyle = 5;
35
36Color_t objColor = kBlack;
37Color_t trackColor = kBlack;
38Color_t towerColor = kBlack;
39
40double effLegXmin = 0.22;
41double effLegXmax = 0.5;
42double effLegYmin = 0.22;
43double effLegYmax = 0.4;
44
45struct resolPlot
46{
47 TH1 *cenResolHist;
48 TH1 *fwdResolHist;
49 double ptmin;
50 double ptmax;
51 double xmin;
52 double xmax;
53 TString obj;
54
55 resolPlot();
56 resolPlot(double ptdown, double ptup, TString object);
57 void set(double ptdown, double ptup, TString object, double xmin = 0, double xmax = 2);
58 void print(){std::cout << ptmin << std::endl;}
59};
60
61
62resolPlot::resolPlot()
63{
64}
65
66resolPlot::resolPlot(double ptdown, double ptup, TString object)
67{
68 this->set(ptdown,ptup,object);
69}
70
71void resolPlot::set(double ptdown, double ptup, TString object, double xmin, double xmax){
72 ptmin = ptdown;
73 ptmax = ptup;
74 obj = object;
75
76 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", 500, xmin, xmax);
77 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", 500, xmin, xmax);
78
79}
80
81void HistogramsCollection(std::vector<resolPlot> *histos, double ptmin, double ptmax, TString obj, double xmin = 0, double xmax = 2)
82{
83 double width;
84 double ptdown;
85 double ptup;
86 resolPlot ptemp;
87
88 for (int i = 0; i < Nbins; i++)
89 {
90 width = (ptmax - ptmin) / Nbins;
91 ptdown = TMath::Power(10,ptmin + i * width );
92 ptup = TMath::Power(10,ptmin + (i+1) * width );
93 ptemp.set(ptdown, ptup, obj, xmin, xmax);
94 histos->push_back(ptemp);
95 }
96}
97
98//------------------------------------------------------------------------------
99
100class ExRootResult;
101class ExRootTreeReader;
102
103//------------------------------------------------------------------------------
104
105void BinLogX(TH1*h)
106{
107
108 TAxis *axis = h->GetXaxis();
109 int bins = axis->GetNbins();
110
111 Axis_t from = axis->GetXmin();
112 Axis_t to = axis->GetXmax();
113 Axis_t width = (to - from) / bins;
114 Axis_t *new_bins = new Axis_t[bins + 1];
115
116 for (int i = 0; i <= bins; i++) {
117 new_bins[i] = TMath::Power(10, from + i * width);
118
119 }
120 axis->Set(bins, new_bins);
121 delete new_bins;
122}
123
124
125//------------------------------------------------------------------------------
126
127template<typename T>
128std::pair<TH1D*, TH1D*> GetEff(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, ExRootTreeReader *treeReader)
129{
130
131 cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
132
133 Long64_t allEntries = treeReader->GetEntries();
134
135 GenParticle *particle;
136 T *recoObj;
137
138 TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
139
140 Float_t deltaR;
141 Float_t pt, eta;
142 Long64_t entry;
143
144 Int_t i, j;
145
146 TH1D *histGenPt = new TH1D(name+" gen spectra Pt",name+" gen spectra Pt", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
147 TH1D *histRecoPt = new TH1D(name+" reco spectra Pt",name+" reco spectra Pt", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
148 TH1D *histGenEta = new TH1D(name+" gen spectra Eta",name+" gen spectra Eta", 12, -3, 3);
149 TH1D *histRecoEta = new TH1D(name+" reco spectra Eta",name+" reco spectra Eta", 12, -3, 3);
150
151
152 BinLogX(histGenPt);
153 BinLogX(histRecoPt);
154
155 // Loop over all events
156 for(entry = 0; entry < allEntries; ++entry)
157 {
158 // Load selected branches with data from specified event
159 treeReader->ReadEntry(entry);
160
161 // Loop over all generated particle in event
162 for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
163 {
164
165 particle = (GenParticle*) branchParticle->At(i);
166 genMomentum = particle->P4();
167
168 deltaR = 999;
169
170 if (particle->PID == pdgID && genMomentum.Pt() > ptrangemin && genMomentum.Pt() < ptrangemax )
171 {
172
173 // Loop over all reco object in event
174 for(j = 0; j < branchReco->GetEntriesFast(); ++j)
175 {
176 recoObj = (T*)branchReco->At(j);
177 recoMomentum = recoObj->P4();
178 // this is simply to avoid warnings from initial state particle
179 // having infite rapidity ...
180 //if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
181
182 // take the closest parton candidate
183 if(TMath::Abs(pdgID) == 5)
184 {
185 Jet *jet = (Jet *)recoObj;
186 if(jet->BTag != 1) continue;
187 }
188 if(TMath::Abs(pdgID) == 15)
189 {
190 Jet *jet = (Jet *)recoObj;
191 if(jet->TauTag != 1) continue;
192 }
193 if(genMomentum.DeltaR(recoMomentum) < deltaR)
194 {
195 deltaR = genMomentum.DeltaR(recoMomentum);
196 bestRecoMomentum = recoMomentum;
197 }
198 }
199
200 pt = genMomentum.Pt();
201 eta = genMomentum.Eta();
202
203 histGenPt->Fill(pt);
204 histGenEta->Fill(eta);
205
206 if(deltaR < 0.3)
207 {
208 histRecoPt->Fill(pt);
209 histRecoEta->Fill(eta);
210 }
211 }
212 }
213 }
214
215
216 std::pair<TH1D*,TH1D*> histos;
217
218 histRecoPt->Divide(histGenPt);
219 histRecoEta->Divide(histGenEta);
220
221 histos.first = histRecoPt;
222 histos.second = histRecoEta;
223
224 return histos;
225}
226
227template<typename T>
228void GetEres(std::vector<resolPlot> *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, ExRootTreeReader *treeReader)
229{
230 Long64_t allEntries = treeReader->GetEntries();
231
232 cout << "** Computing resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
233
234 GenParticle *particle;
235 T* recoObj;
236
237 TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
238
239 Float_t deltaR;
240 Float_t pt, eta;
241 Long64_t entry;
242
243 Int_t i, j, bin;
244
245 // Loop over all events
246 for(entry = 0; entry < allEntries; ++entry)
247 {
248 // Load selected branches with data from specified event
249 treeReader->ReadEntry(entry);
250
251 // Loop over all reconstructed jets in event
252 for(i = 0; i < branchReco->GetEntriesFast(); ++i)
253 {
254 recoObj = (T*) branchReco->At(i);
255 recoMomentum = recoObj->P4();
256
257 deltaR = 999;
258
259 // Loop over all hard partons in event
260 for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
261 {
262 particle = (GenParticle*) branchParticle->At(j);
263 if (particle->PID == pdgID && particle->Status == 1)
264 {
265 genMomentum = particle->P4();
266
267 // this is simply to avoid warnings from initial state particle
268 // having infite rapidity ...
269 if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
270
271 // take the closest parton candidate
272 if(genMomentum.DeltaR(recoMomentum) < deltaR)
273 {
274 deltaR = genMomentum.DeltaR(recoMomentum);
275 bestGenMomentum = genMomentum;
276 }
277 }
278 }
279
280 if(deltaR < 0.3)
281 {
282 pt = bestGenMomentum.Pt();
283 eta = TMath::Abs(bestGenMomentum.Eta());
284
285 for (bin = 0; bin < Nbins; bin++)
286 {
287 if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta > 0.0 && eta < 2.5)
288 {
289 if (eta < 1.5) {histos->at(bin).cenResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
290 else if (eta < 2.5) {histos->at(bin).fwdResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
291 }
292 }
293 }
294 }
295 }
296}
297void GetJetsEres(std::vector<resolPlot> *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader)
298{
299
300 Long64_t allEntries = treeReader->GetEntries();
301
302 cout << "** Computing resolution of " << branchJet->GetName() << " induced by " << branchGenJet->GetName() << endl;
303
304 Jet *jet, *genjet;
305
306 TLorentzVector jetMomentum, genJetMomentum, bestGenJetMomentum;
307
308 Float_t deltaR;
309 Float_t pt, eta;
310 Long64_t entry;
311
312 Int_t i, j, bin;
313
314 // Loop over all events
315 for(entry = 0; entry < allEntries; ++entry)
316 {
317 // Load selected branches with data from specified event
318 treeReader->ReadEntry(entry);
319
320 if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
321
322 // Loop over all reconstructed jets in event
323 for(i = 0; i < TMath::Min(2,branchJet->GetEntriesFast()); ++i) //branchJet->GetEntriesFast(); ++i)
324 {
325
326 jet = (Jet*) branchJet->At(i);
327 jetMomentum = jet->P4();
328
329 deltaR = 999;
330
331 // Loop over all hard partons in event
332 for(j = 0; j < TMath::Min(2,branchGenJet->GetEntriesFast()); ++j)
333 {
334 genjet = (Jet*) branchGenJet->At(j);
335
336 genJetMomentum = genjet->P4();
337
338 // this is simply to avoid warnings from initial state particle
339 // having infite rapidity ...
340 if(genJetMomentum.Px() == 0 && genJetMomentum.Py() == 0) continue;
341
342 // take the closest parton candidate
343 if(genJetMomentum.DeltaR(jetMomentum) < deltaR)
344 {
345 deltaR = genJetMomentum.DeltaR(jetMomentum);
346 bestGenJetMomentum = genJetMomentum;
347 }
348 }
349
350 if(deltaR < 0.25)
351 {
352 pt = genJetMomentum.Pt();
353 eta = TMath::Abs(genJetMomentum.Eta());
354
355 for (bin = 0; bin < Nbins; bin++)
356 {
357 if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta > 0.0 && eta < 1.5)
358 {
359 histos->at(bin).cenResolHist->Fill(jetMomentum.Pt()/bestGenJetMomentum.Pt());
360 }
361 }
362 }
363 }
364 }
365}
366
367void GetMetres(std::vector<resolPlot> *histos, TClonesArray *branchScalarHT, TClonesArray *branchMet, TClonesArray *branchJet, ExRootTreeReader *treeReader)
368{
369
370 Long64_t allEntries = treeReader->GetEntries();
371
372 cout << "** Computing resolution of " << branchMet->GetName() << " vs " << branchScalarHT->GetName() << endl;
373
374 MissingET *met;
375 ScalarHT *scalarHT;
376
377 Long64_t entry;
378
379 Int_t bin;
380 Double_t ht;
381
382 // Loop over all events
383 for(entry = 0; entry < allEntries; ++entry)
384 {
385 // Load selected branches with data from specified event
386 treeReader->ReadEntry(entry);
387
388 if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
389
390 if (branchJet->GetEntriesFast() > 1)
391 {
392 met = (MissingET*) branchMet->At(0);
393 scalarHT = (ScalarHT*) branchScalarHT->At(0);
394 ht = scalarHT->HT;
395
396 for (bin = 0; bin < Nbins; bin++)
397 {
398 if(ht > histos->at(bin).ptmin && ht < histos->at(bin).ptmax )
399 {
400 histos->at(bin).cenResolHist->Fill(met->P4().Px());
401 histos->at(bin).fwdResolHist->Fill(met->P4().Py());
402 }
403 }
404 }
405 }
406}
407
408
409std::pair<Double_t, Double_t> GausFit(TH1* hist)
410{
411 TF1 *f1 = new TF1("f1", "gaus", hist->GetMean()-2*hist->GetRMS(), hist->GetMean()+2*hist->GetRMS());
412 hist->Fit("f1","RQ");
413 Double_t sig = f1->GetParameter(2);
414 Double_t sigErr = f1->GetParError(2);
415 delete f1;
416 return make_pair (sig, sigErr);
417 //return make_pair (hist->GetRMS(), hist->GetRMSError());
418}
419
420
421TGraphErrors EresGraph(std::vector<resolPlot> *histos, bool central, bool rms = false)
422{
423 Int_t bin;
424 Int_t count = 0;
425 TGraphErrors gr = TGraphErrors(Nbins/2);
426 Double_t sig = 0;
427 Double_t sigErr = 0;
428 for (bin = 0; bin < Nbins; bin++)
429 {
430 if (central == true && histos->at(bin).cenResolHist->GetEntries() > 100)
431 {
432 std::cout << " pt : " << (histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0;
433 std::cout << " mean : " << histos->at(bin).cenResolHist->GetMean() << " RMS : " << histos->at(bin).cenResolHist->GetRMS();
434 std::cout << " entries : " << histos->at(bin).cenResolHist->GetEntries() << std::endl;
435 std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).cenResolHist);
436 if (rms == true)
437 {
438 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
439 gr.SetPointError(count,0, sigvalues.second); // to correct
440 }
441 else
442 {
443 gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
444 gr.SetPointError(count,0, sigvalues.second);
445 }
446 count++;
447 }
448 /*
449 else if (histos->at(bin).cenResolHist->GetEntries() > 10)
450 {
451 histos->at(bin).fwdResolHist->Fit("gaus","","", -2*histos->at(bin).fwdResolHist->GetRMS(), 2*histos->at(bin).fwdResolHist->GetRMS());
452 TF1 *f = histos->at(bin).fwdResolHist->GetFunction("gaus");
453 Double_t sig = f->GetParameter(2);
454 Double_t sigErr = f->GetParError(2);
455 gr.SetPoint(bin,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sig);
456 gr.SetPointError(bin,0, sigErr);
457 }
458 */
459 }
460 return gr;
461}
462
463
464//------------------------------------------------------------------------------
465
466
467// type 1 : object, 2 : track, 3 : tower
468
469void addGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int type)
470{
471
472 gr->SetLineWidth(3);
473
474 switch ( type )
475 {
476 case 1:
477 gr->SetLineColor(objColor);
478 gr->SetLineStyle(objStyle);
479 std::cout << "Adding " << gr->GetName() << std::endl;
480 mg->Add(gr);
481 leg->AddEntry(gr,"Reco","l");
482 break;
483
484 case 2:
485 gr->SetLineColor(trackColor);
486 gr->SetLineStyle(trackStyle);
487 mg->Add(gr);
488 leg->AddEntry(gr,"Track","l");
489 break;
490
491 case 3:
492 gr->SetLineColor(towerColor);
493 gr->SetLineStyle(towerStyle);
494 mg->Add(gr);
495 leg->AddEntry(gr,"Tower","l");
496 break;
497
498 default:
499 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
500 break;
501 }
502}
503
504void addHist(TH1D *h, TLegend *leg, int type)
505{
506 h->SetLineWidth(3);
507
508 switch ( type )
509 {
510 case 1:
511 h->SetLineColor(objColor);
512 h->SetLineStyle(objStyle);
513 leg->AddEntry(h,"Reco","l");
514 break;
515
516 case 2:
517 h->SetLineColor(trackColor);
518 h->SetLineStyle(trackStyle);
519 leg->AddEntry(h,"Track","l");
520 break;
521
522 case 3:
523 h->SetLineColor(towerColor);
524 h->SetLineStyle(towerStyle);
525 leg->AddEntry(h,"Tower","l");
526 break;
527
528 default:
529 std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
530 break;
531 }
532}
533
534void DrawAxis(TMultiGraph *mg, TLegend *leg, double max)
535{
536 mg->SetMinimum(0.);
537 mg->SetMaximum(max);
538 mg->GetXaxis()->SetLimits(ptrangemin,ptrangemax);
539 mg->GetYaxis()->SetTitle("#sigma (E) / E_{gen}");
540 mg->GetXaxis()->SetTitle("E_{gen}");
541 mg->GetYaxis()->SetTitleSize(0.07);
542 mg->GetXaxis()->SetTitleSize(0.07);
543 mg->GetYaxis()->SetLabelSize(0.07);
544 mg->GetXaxis()->SetLabelSize(0.07);
545
546 leg->SetLineStyle(0);
547 leg->SetFillStyle(0);
548 leg->SetLineWidth(0);
549 leg->SetLineColor(0);
550
551 gStyle->SetOptTitle(0);
552 gPad->SetLogx();
553 gPad->SetBottomMargin(0.2);
554 gPad->SetLeftMargin(0.2);
555 gPad->Modified();
556 gPad->Update();
557
558}
559
560void DrawAxis(TH1D *h, TLegend *leg, int type)
561{
562
563 h->GetYaxis()->SetRangeUser(0,1.2);
564 if (type == 0) h->GetXaxis()->SetTitle("E_{gen}");
565 else h->GetXaxis()->SetTitle("#eta");
566 h->GetYaxis()->SetTitle("#epsilon");
567 h->GetYaxis()->SetTitleSize(0.07);
568 h->GetXaxis()->SetTitleSize(0.07);
569 h->GetYaxis()->SetLabelSize(0.07);
570 h->GetXaxis()->SetLabelSize(0.07);
571 leg->SetLineStyle(0);
572 leg->SetFillStyle(0);
573 leg->SetLineWidth(0);
574 leg->SetLineColor(0);
575
576 gStyle->SetOptTitle(0);
577 gStyle->SetOptStat(0);
578 gPad->SetBottomMargin(0.2);
579 gPad->SetLeftMargin(0.2);
580
581 gPad->Modified();
582 gPad->Update();
583
584}
585
586
587void Validation(const char *inputFile, const char *outputFile)
588{
589 //gSystem->Load("libDelphes");
590
591 std::cout << "input file : " << inputFile << " " << " , output file : " << outputFile << std::endl;
592
593 TChain *chain = new TChain("Delphes");
594 chain->Add(inputFile);
595
596 ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
597 TClonesArray *branchParticle = treeReader->UseBranch("Particle");
598 TClonesArray *branchElectron = treeReader->UseBranch("Electron");
599 TClonesArray *branchMuon = treeReader->UseBranch("Muon");
600 TClonesArray *branchPhoton = treeReader->UseBranch("Photon");
601 TClonesArray *branchTrack = treeReader->UseBranch("Track");
602 TClonesArray *branchTower = treeReader->UseBranch("Tower");
603 TClonesArray *branchGenJet = treeReader->UseBranch("GenJet");
604 TClonesArray *branchPFJet = treeReader->UseBranch("Jet");
605 TClonesArray *branchCaloJet = treeReader->UseBranch("CaloJet");
606 TClonesArray *branchScalarHT = treeReader->UseBranch("ScalarHT");
607 TClonesArray *branchMet = treeReader->UseBranch("MissingET");
608
609#ifdef ELECTRON
610
611 ///////////////
612 // Electrons //
613 ///////////////
614
615 // Reconstruction efficiency
616 TString elecs = "Electron";
617 int elID = 11;
618 std::pair<TH1D*,TH1D*> histos_el = GetEff<Electron>(branchElectron, branchParticle, "Electron", elID, treeReader);
619
620 // tracking reconstruction efficiency
621 std::pair <TH1D*,TH1D*> histos_eltrack = GetEff<Track>(branchTrack, branchParticle, "electronTrack", elID, treeReader);
622
623 // Tower reconstruction efficiency
624 std::pair <TH1D*,TH1D*> histos_eltower = GetEff<Tower>(branchTower, branchParticle, "electronTower", elID, treeReader);
625
626 // Electron Energy Resolution
627 std::vector<resolPlot> plots_el;
628 HistogramsCollection(&plots_el, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electrons");
629 GetEres<Electron>( &plots_el, branchElectron, branchParticle, elID, treeReader);
630 TGraphErrors gr_el = EresGraph(&plots_el, true);
631 gr_el.SetName("Electron");
632
633 // Electron Track Energy Resolution
634 std::vector<resolPlot> plots_eltrack;
635 HistogramsCollection(&plots_eltrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTracks");
636 GetEres<Track>( &plots_eltrack, branchTrack, branchParticle, elID, treeReader);
637 TGraphErrors gr_eltrack = EresGraph(&plots_eltrack, true);
638 gr_eltrack.SetName("ElectronTracks");
639
640 // Electron Tower Energy Resolution
641 std::vector<resolPlot> plots_eltower;
642 HistogramsCollection(&plots_eltower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTower");
643 GetEres<Tower>( &plots_eltower, branchTower, branchParticle, elID, treeReader);
644 TGraphErrors gr_eltower = EresGraph(&plots_eltower, true);
645 gr_eltower.SetName("ElectronTower");
646
647 // Canvases
648 TString elEff = "electronEff";
649 TCanvas *C_el1 = new TCanvas(elEff,elEff, 1000, 500);
650 C_el1->Divide(2);
651 C_el1->cd(1);
652 TLegend *leg_el1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
653
654 gPad->SetLogx();
655 histos_eltrack.first->Draw();
656 addHist(histos_eltrack.first, leg_el1, 2);
657 histos_eltower.first->Draw("same");
658 addHist(histos_eltower.first, leg_el1, 3);
659 histos_el.first->Draw("same");
660 addHist(histos_el.first, leg_el1, 1);
661
662 DrawAxis(histos_eltrack.first, leg_el1, 0);
663 leg_el1->Draw();
664
665 C_el1->cd(2);
666 TLegend *leg_el2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
667
668 histos_eltrack.second->Draw();
669 addHist(histos_eltrack.second, leg_el2, 2);
670 histos_eltower.second->Draw("same");
671 addHist(histos_eltower.second, leg_el2, 3);
672 histos_el.second->Draw("same");
673 addHist(histos_el.second, leg_el2, 1);
674
675 DrawAxis(histos_eltrack.second, leg_el2, 1);
676 delete(leg_el2);
677 C_el1->cd(0);
678
679 TString elRes = "electronERes";
680 TCanvas *C_el2 = new TCanvas(elRes,elRes, 1000, 500);
681 TMultiGraph *mg_el = new TMultiGraph(elRes,elRes);
682 TLegend *leg_el = new TLegend(0.52,0.7,0.9,0.9);
683
684 addGraph(mg_el, &gr_eltower, leg_el, 3);
685 addGraph(mg_el, &gr_eltrack, leg_el, 2);
686 addGraph(mg_el, &gr_el, leg_el, 1);
687
688 mg_el->Draw("ACX");
689 leg_el->Draw();
690
691 DrawAxis(mg_el, leg_el, 0.1);
692
693 C_el1->SaveAs(elEff+".eps");
694 C_el2->SaveAs(elRes+".eps");
695
696#endif
697
698 gDirectory->cd(0);
699
700#ifdef MUON
701
702 ///////////
703 // Muons //
704 ///////////
705
706 // Reconstruction efficiency
707 int muID = 13;
708 std::pair<TH1D*,TH1D*> histos_mu = GetEff<Muon>(branchMuon, branchParticle,"Muon", muID, treeReader);
709
710 // muon tracking reconstruction efficiency
711 std::pair <TH1D*,TH1D*> histos_mutrack = GetEff<Track>(branchTrack, branchParticle, "muonTrack", muID, treeReader);
712
713 // Muon Energy Resolution
714 std::vector<resolPlot> plots_mu;
715 HistogramsCollection(&plots_mu, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muons");
716 GetEres<Muon>( &plots_mu, branchMuon, branchParticle, muID, treeReader);
717 TGraphErrors gr_mu = EresGraph(&plots_mu, true);
718 gr_mu.SetName("Muon");
719
720 // Muon Track Energy Resolution
721 std::vector<resolPlot> plots_mutrack;
722 HistogramsCollection(&plots_mutrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muonsTracks");
723 GetEres<Track>( &plots_mutrack, branchTrack, branchParticle, muID, treeReader);
724 TGraphErrors gr_mutrack = EresGraph(&plots_mutrack, true);
725 gr_mutrack.SetName("MuonTracks");
726
727 // Canvas
728
729 TString muEff = "muonEff";
730 TCanvas *C_mu1 = new TCanvas(muEff,muEff, 1000, 500);
731 C_mu1->Divide(2);
732 C_mu1->cd(1);
733 TLegend *leg_mu1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
734
735 gPad->SetLogx();
736 histos_mutrack.first->Draw();
737 addHist(histos_mutrack.first, leg_mu1, 2);
738 histos_mu.first->Draw("same");
739 addHist(histos_mu.first, leg_mu1, 1);
740
741 DrawAxis(histos_mutrack.first, leg_mu1, 0);
742 leg_mu1->Draw();
743
744 C_mu1->cd(2);
745 TLegend *leg_mu2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
746
747 histos_mutrack.second->Draw();
748 addHist(histos_mutrack.second, leg_mu2, 2);
749 histos_mu.second->Draw("same");
750 addHist(histos_mu.second, leg_mu2, 1);
751
752 DrawAxis(histos_mutrack.second, leg_mu2, 1);
753
754 TString muRes = "muonERes";
755 TCanvas *C_mu = new TCanvas(muRes,muRes, 1000, 500);
756 TMultiGraph *mg_mu = new TMultiGraph(muRes,muRes);
757 TLegend *leg_mu = new TLegend(0.52,0.7,0.9,0.9);
758
759 addGraph(mg_mu, &gr_mutrack, leg_mu, 2);
760 addGraph(mg_mu, &gr_mu, leg_mu, 1);
761
762 mg_mu->Draw("ACX");
763 leg_mu->Draw();
764
765 DrawAxis(mg_mu, leg_mu, 0.3);
766
767 C_mu1->SaveAs(muEff+".eps");
768 C_mu->SaveAs(muRes+".eps");
769
770#endif
771
772 gDirectory->cd(0);
773
774#ifdef PHOTON
775
776 /////////////
777 // Photons //
778 /////////////
779
780 // Reconstruction efficiency
781 int phID = 22;
782 std::pair<TH1D*,TH1D*> histos_ph = GetEff<Electron>(branchPhoton, branchParticle, "Photon", phID, treeReader);
783 std::pair<TH1D*,TH1D*> histos_phtower = GetEff<Electron>(branchTower, branchParticle, "Photon", phID, treeReader);
784
785 // Photon Energy Resolution
786 std::vector<resolPlot> plots_ph;
787 HistogramsCollection(&plots_ph, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photons");
788 GetEres<Photon>( &plots_ph, branchPhoton, branchParticle, phID, treeReader);
789 TGraphErrors gr_ph = EresGraph(&plots_ph, true);
790 gr_ph.SetName("Photon");
791
792 // Photon Tower Energy Resolution
793 std::vector<resolPlot> plots_phtower;
794 HistogramsCollection(&plots_phtower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photonsTower");
795 GetEres<Tower>( &plots_phtower, branchTower, branchParticle, phID, treeReader);
796 TGraphErrors gr_phtower = EresGraph(&plots_phtower, true);
797 gr_phtower.SetName("PhotonTower");
798
799 // Canvas
800
801 TString phEff = "photonEff";
802 TCanvas *C_ph1 = new TCanvas(phEff,phEff, 1000, 500);
803 C_ph1->Divide(2);
804 C_ph1->cd(1);
805 TLegend *leg_ph1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
806
807 gPad->SetLogx();
808 histos_phtower.first->Draw();
809 addHist(histos_phtower.first, leg_ph1, 3);
810 histos_ph.first->Draw("same");
811 addHist(histos_ph.first, leg_ph1, 1);
812
813 DrawAxis(histos_phtower.first, leg_ph1, 0);
814 leg_ph1->Draw();
815
816 C_ph1->cd(2);
817 TLegend *leg_ph2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
818
819 histos_phtower.second->Draw("same");
820 addHist(histos_phtower.second, leg_ph2, 3);
821 histos_ph.second->Draw("same");
822 addHist(histos_ph.second, leg_ph2, 1);
823
824 DrawAxis(histos_phtower.second, leg_ph2, 1);
825
826 C_ph1->SaveAs(phEff+".eps");
827
828 TString phRes = "phERes";
829 TCanvas *C_ph = new TCanvas(phRes,phRes, 1000, 500);
830 TMultiGraph *mg_ph = new TMultiGraph(phRes,phRes);
831 TLegend *leg_ph = new TLegend(0.52,0.7,0.9,0.9);
832
833 addGraph(mg_ph, &gr_phtower, leg_ph, 3);
834 addGraph(mg_ph, &gr_ph, leg_ph, 1);
835
836 mg_ph->Draw("ACX");
837 leg_ph->Draw();
838
839 DrawAxis(mg_ph, leg_ph, 0.3);
840
841 C_ph->SaveAs(phRes+".eps");
842
843#endif
844
845 gDirectory->cd(0);
846
847#ifdef JET
848
849 //////////
850 // Jets //
851 //////////
852
853 // BJets Reconstruction efficiency
854 int bID = 5;
855 std::pair<TH1D*,TH1D*> histos_btag = GetEff<Jet>(branchPFJet, branchParticle,"BTag", bID, treeReader);
856
857 // TauJets Reconstruction efficiency
858 int tauID = 15;
859 std::pair<TH1D*,TH1D*> histos_tautag = GetEff<Jet>(branchPFJet, branchParticle,"TauTag", tauID, treeReader);
860
861 // PFJets Energy Resolution
862 std::vector<resolPlot> plots_pfjets;
863 HistogramsCollection(&plots_pfjets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "PFJet");
864 GetJetsEres( &plots_pfjets, branchPFJet, branchGenJet, treeReader);
865 TGraphErrors gr_pfjets = EresGraph(&plots_pfjets, true);
866 gr_pfjets.SetName("pfJet");
867
868 // CaloJets Energy Resolution
869 std::vector<resolPlot> plots_calojets;
870 HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "CaloJet");
871 GetJetsEres( &plots_calojets, branchCaloJet, branchGenJet, treeReader);
872 TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
873 gr_calojets.SetName("caloJet");
874
875 // MET Resolution vs HT
876 std::vector<resolPlot> plots_met;
877 HistogramsCollection(&plots_met, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "MET", -200, 200);
878 GetMetres( &plots_met, branchScalarHT, branchMet, branchPFJet, treeReader);
879 TGraphErrors gr_met = EresGraph(&plots_met, true);
880 gr_calojets.SetName("MET");
881
882 // Canvas
883
884 TString btagEff = "btagEff";
885 TCanvas *C_btag1 = new TCanvas(btagEff,btagEff, 1000, 500);
886 C_btag1->Divide(2);
887 C_btag1->cd(1);
888 TLegend *leg_btag1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
889
890 gPad->SetLogx();
891 histos_btag.first->Draw();
892 addHist(histos_btag.first, leg_btag1, 1);
893
894 DrawAxis(histos_btag.first, leg_btag1, 0);
895 leg_btag1->Draw();
896
897 C_btag1->cd(2);
898 TLegend *leg_btag2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
899
900 histos_btag.second->Draw();
901 addHist(histos_btag.second, leg_btag2, 1);
902
903 DrawAxis(histos_btag.second, leg_btag2, 1);
904 delete(leg_btag2);
905 C_btag1->cd(0);
906
907 TString tautagEff = "tautagEff";
908 TCanvas *C_tautag1 = new TCanvas(tautagEff,tautagEff, 1000, 500);
909 C_tautag1->Divide(2);
910 C_tautag1->cd(1);
911 TLegend *leg_tautag1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
912
913 gPad->SetLogx();
914 histos_tautag.first->Draw();
915 addHist(histos_tautag.first, leg_tautag1, 1);
916
917 DrawAxis(histos_tautag.first, leg_tautag1, 0);
918 leg_tautag1->Draw();
919
920 C_tautag1->cd(2);
921 TLegend *leg_tautag2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
922
923 histos_tautag.second->Draw();
924 addHist(histos_tautag.second, leg_tautag2, 1);
925
926 DrawAxis(histos_tautag.second, leg_tautag2, 1);
927 delete(leg_tautag2);
928 C_tautag1->cd(0);
929
930 TString jetRes = "jetERes";
931 TCanvas *C_jet = new TCanvas(jetRes,jetRes, 1000, 500);
932 TMultiGraph *mg_jet = new TMultiGraph(jetRes,jetRes);
933 TLegend *leg_jet = new TLegend(0.52,0.7,0.9,0.9);
934
935 addGraph(mg_jet, &gr_calojets, leg_jet, 3);
936 addGraph(mg_jet, &gr_pfjets, leg_jet, 1);
937
938 mg_jet->Draw("ACX");
939 leg_jet->Draw();
940
941 DrawAxis(mg_jet, leg_jet, 0.25);
942
943 C_btag1->SaveAs(btagEff+".eps");
944 C_jet->SaveAs(jetRes+".eps");
945
946 TString metRes = "MetRes";
947 TCanvas *C_met = new TCanvas(metRes,metRes, 1000, 500);
948
949 TMultiGraph *mg_met = new TMultiGraph(metRes,metRes);
950 TLegend *leg_met = new TLegend(0.52,0.7,0.9,0.9);
951
952 addGraph(mg_met, &gr_met, leg_met, 3);
953
954 mg_met->Draw("ACX");
955 leg_met->Draw();
956
957 DrawAxis(mg_met, leg_met, 100);
958
959 C_met->SaveAs(metRes+".eps");
960
961
962#endif
963
964 /*
965 // CaloJets Energy Resolution
966 std::vector<resolPlot> plots_calojets;
967 HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "caloJet");
968 GetJetsEres( &plots_calojets, branchCaloJet, branchGenJet, treeReader);
969 TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
970 gr_calojets.SetName("caloJet");
971 */
972
973 TFile *fout = new TFile(outputFile,"recreate");
974
975#ifdef ELECTRON
976 for (int bin = 0; bin < Nbins; bin++)
977 {
978 plots_el.at(bin).cenResolHist->Write();
979 plots_eltrack.at(bin).cenResolHist->Write();
980 plots_eltower.at(bin).cenResolHist->Write();
981 }
982 // gr.Write();
983 histos_el.first->Write();
984 //histos_el.second->Write();
985 histos_eltrack.first->Write();
986 //histos_eltrack.second->Write();
987 histos_eltower.first->Write();
988 C_el1->Write();
989 C_el2->Write();
990#endif
991
992#ifdef MUON
993 histos_mu.first->Write();
994 histos_mu.second->Write();
995 histos_mutrack.first->Write();
996 histos_mutrack.second->Write();
997 C_mu1->Write();
998 C_mu->Write();
999#endif
1000
1001#ifdef PHOTON
1002 histos_ph.first->Write();
1003 histos_ph.second->Write();
1004 C_ph1->Write();
1005 C_ph->Write();
1006#endif
1007
1008#ifdef JET
1009 for (int bin = 0; bin < Nbins; bin++)
1010 {
1011 plots_pfjets.at(bin).cenResolHist->Write();
1012 plots_calojets.at(bin).cenResolHist->Write();
1013 plots_met.at(bin).cenResolHist->Write();
1014 }
1015 histos_btag.first->Write();
1016 histos_btag.second->Write();
1017 histos_tautag.first->Write();
1018 histos_tautag.second->Write();
1019 C_btag1->Write();
1020 C_tautag1->Write();
1021 C_jet->Write();
1022 C_met->Write();
1023#endif
1024
1025 fout->Write();
1026
1027 cout << "** Exiting..." << endl;
1028
1029 //delete plots;
1030 //delete result;
1031 delete treeReader;
1032 delete chain;
1033}
1034
1035//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.