Fork me on GitHub

source: git/examples/Validation.C@ 734b267

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

validation progress, jets , e, gamma, mu

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