Fork me on GitHub

source: svn/trunk/Resolutions_ATLAS.cpp@ 468

Last change on this file since 468 was 467, checked in by Xavier Rouby, 15 years ago

cleaning, comments, progressbar, cout

File size: 18.5 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \------------- arXiv:0903.2225v1 ------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** ROOT: Nucl. Inst. & Meth. in Phys. Res. A389 (1997) 81-86 **
12** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
13** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
14** FROG: [hep-ex/0901.2718v1] **
15** HepMC: Comput. Phys. Commun.134 (2001) 41 **
16** **
17** ------------------------------------------------------------------ **
18** **
19** Main authors: **
20** ------------- **
21** **
22** Severine Ovyn Xavier Rouby **
23** severine.ovyn@uclouvain.be xavier.rouby@cern **
24** **
25** Center for Particle Physics and Phenomenology (CP3) **
26** Universite catholique de Louvain (UCL) **
27** Louvain-la-Neuve, Belgium **
28** **
29** Copyright (C) 2008-2009, **
30** All rights reserved. **
31** **
32***********************************************************************/
33
34/// \file Resolution_ATLAS.cpp
35/// \brief Resolution for Atlas
36
37#include "TChain.h"
38#include "TApplication.h"
39#include "TFile.h"
40#include "TStopwatch.h"
41
42#include "ExRootTreeReader.h"
43#include "ExRootTreeWriter.h"
44#include "ExRootTreeBranch.h"
45#include "ExRootProgressBar.h"
46#include "TreeClasses.h"
47
48#include "SmearUtil.h"
49#include "JetsUtil.h"
50#include "BFieldProp.h"
51
52#include<vector>
53#include<iostream>
54
55using namespace std;
56
57//------------------------------------------------------------------------------
58
59// //********************************** PYTHIA INFORMATION*********************************
60
61TSimpleArray<TRootGenParticle> TauHadr(const TClonesArray *GEN)
62 {
63 TIter it((TCollection*)GEN);
64 it.Reset();
65 TRootGenParticle *gen1;
66 TSimpleArray<TRootGenParticle> array,array2;
67
68 while((gen1 = (TRootGenParticle*) it.Next()))
69 {
70 array.Add(gen1);
71 }
72 it.Reset();
73 bool tauhad;
74 while((gen1 = (TRootGenParticle*) it.Next()))
75 {
76 tauhad=false;
77 if(abs(gen1->PID)==15)
78 {
79 int d1=gen1->D1;
80 int d2=gen1->D2;
81 if((d1 < array.GetEntries()) && (d1 > 0) && (d2 < array.GetEntries()) && (d2 > 0))
82 {
83 tauhad=true;
84 for(int d=d1; d < d2+1; d++)
85 {
86 if(abs(array[d]->PID)== pE || abs(array[d]->PID)== pMU)tauhad=false;
87 }
88 }
89 }
90 if(tauhad)array2.Add(gen1);
91 }
92 return array2;
93 }
94
95double EnergySmallCone(const vector<TLorentzVector> &towers, const float eta, const float phi,float energy_scone,float JET_seed) {
96 double Energie=0;
97 for(unsigned int i=0; i < towers.size(); i++) {
98 if(towers[i].Pt() < JET_seed) continue;
99 if((DeltaR(phi,eta,towers[i].Phi(),towers[i].Eta()) < energy_scone)) {
100 Energie += towers[i].E();
101 }
102 }
103 return Energie;
104}
105
106
107void PairingJet(TLorentzVector &JETSm, const TLorentzVector &JET, const TClonesArray *branchJet, const float dR=0.25)
108{
109 JETSm.SetPxPyPzE(0,0,0,0);
110 float deltaRtest=5000;
111 TIter itJet((TCollection*)branchJet);
112 TRootJet *jet;
113 itJet.Reset();
114 while( (jet = (TRootJet*) itJet.Next()) )
115 {
116 TLorentzVector Att;
117 Att.SetPtEtaPhiE(jet->PT,jet->Eta,jet->Phi,jet->E);
118 if(DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta()) < deltaRtest)
119 {
120 deltaRtest = DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta());
121 if(deltaRtest < dR)
122 {
123 JETSm = Att;
124 }
125 }
126 }
127}
128
129void PairingElec(TLorentzVector &ELECSm, const TLorentzVector &ELEC, const TClonesArray *branchElec)
130{
131 ELECSm.SetPxPyPzE(0,0,0,0);
132 float deltaRtest=5000;
133 TIter itElec((TCollection*)branchElec);
134 TRootElectron *elec;
135 itElec.Reset();
136 while( (elec = (TRootElectron*) itElec.Next()) )
137 {
138 TLorentzVector Att;
139 Att.SetPtEtaPhiE(elec->PT,elec->Eta,elec->Phi,elec->E);
140 if(DeltaR(ELEC.Phi(),ELEC.Eta(),Att.Phi(),Att.Eta()) < deltaRtest)
141 {
142 deltaRtest = DeltaR(ELEC.Phi(),ELEC.Eta(),Att.Phi(),Att.Eta());
143 if(deltaRtest < 0.025)
144 {
145 ELECSm = Att;
146 }
147 }
148 }
149}
150
151void PairingMuon(TLorentzVector &MUONSm, const TLorentzVector &MUON, const TClonesArray *branchMuon)
152{
153 MUONSm.SetPxPyPzE(0,0,0,0);
154 float deltaRtest=5000;
155 TIter itMuon((TCollection*)branchMuon);
156 TRootMuon *muon;
157 itMuon.Reset();
158 while( (muon = (TRootMuon*) itMuon.Next()) )
159 {
160 TLorentzVector Att;
161 Att.SetPxPyPzE(muon->Px,muon->Py,muon->Pz,muon->E);
162 if(DeltaR(MUON.Phi(),MUON.Eta(),Att.Phi(),Att.Eta()) < deltaRtest)
163 {
164 deltaRtest = DeltaR(MUON.Phi(),MUON.Eta(),Att.Phi(),Att.Eta());
165 if(deltaRtest < 0.025)
166 {
167 MUONSm = Att;
168 }
169 }
170 }
171}
172
173unsigned int NumTracks(const TClonesArray *branchTracks, const float pt_track, const float eta, const float phi,float track_scone) {
174 unsigned int numtrack=0;
175 TIter itTrack((TCollection*)branchTracks);
176 TRootTracks *track;
177 itTrack.Reset();
178 while( (track = (TRootTracks*) itTrack.Next()) )
179 {
180 if((track->PT < pt_track )||
181 (DeltaR(phi,eta,track->Phi,track->Eta) > track_scone)
182 )continue;
183 numtrack++;
184 }
185 return numtrack;
186}
187
188
189
190int main(int argc, char *argv[])
191{
192 int appargc = 2;
193 char *appName= new char[20];
194 char *appOpt= new char[20];
195 sprintf(appName,"Resolution_ATLAS");
196 sprintf(appOpt,"-b");
197 char *appargv[] = {appName,appOpt};
198 TApplication app(appName, &appargc, appargv);
199 delete [] appName;
200 delete [] appOpt;
201
202 if(argc != 3) {
203 cout << " Usage: " << argv[0] << " input_file" << " output_file" << endl;
204 cout << " input_file - input file in Delphes-root format," << endl;
205 cout << " output_file - output file." << endl;
206 exit(1);
207 }
208
209
210 // 1. ********** initialisation ***********
211
212 srand (time (NULL));
213 TStopwatch globalwatch;
214 globalwatch.Start();
215
216 string inputfilename(argv[1]), outputfilename(argv[2]);
217
218 // 1.1 checks the input file
219 if(inputfilename.find(".root") > inputfilename.length() ) {
220 cout << "input_file should be a .root file from Delphes!\n";
221 return -1;
222 }
223 TFile f(inputfilename.c_str());
224 if (!f.FindKey("GEN") || !f.FindKey("Analysis") ) {
225 cout << "input_file should be a .root file from Delphes!\n";
226 cout << "it should contain both \"GEN\" and \"Analysis\" trees at least.\n";
227 return -1;
228 }
229 f.Close();
230
231 // 1.2 checks the output file
232 if(outputfilename.find(".root") > outputfilename.length() ) {
233 cout << "output_file should be a .root file!\n";
234 return -1;
235 }
236 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE");// Creates the file, but should be closed just after
237 outputFile->Close();
238
239 // 1.3 Reads the trees in input file
240 TChain chainGEN("GEN");
241 chainGEN.Add(inputfilename.c_str());
242 ExRootTreeReader *treeReaderGEN = new ExRootTreeReader(&chainGEN);
243 TChain chain("Analysis");
244 chain.Add(inputfilename.c_str());
245 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
246 const TClonesArray *branchJet = treeReader->UseBranch("Jet");
247 const TClonesArray *branchElec = treeReader->UseBranch("Electron");
248 const TClonesArray *branchMuon = treeReader->UseBranch("Muon");
249 const TClonesArray *branchTracks = treeReader->UseBranch("Tracks");
250 const TClonesArray *branchTowers = treeReader->UseBranch("CaloTower");
251 const TClonesArray *branchGen = treeReaderGEN->UseBranch("Particle");
252 TIter itGen((TCollection*)branchGen);
253
254
255 // 1.4 Prepares the output root file
256 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
257 ExRootTreeBranch *branchjet = treeWriter->NewBranch("JetPTResol", RESOLJET::Class());
258 ExRootTreeBranch *branchelec = treeWriter->NewBranch("ElecEResol", RESOLELEC::Class());
259 ExRootTreeBranch *branchmuon = treeWriter->NewBranch("MuonPTResol", RESOLMUON::Class());
260 ExRootTreeBranch *branchtaujet = treeWriter->NewBranch("TauJetPTResol", TAUHAD::Class());
261 ExRootTreeBranch *branchetmis = treeWriter->NewBranch("ETmisResol",ETMIS::Class());
262
263 // 1.5 other initialisations
264 TRootGenParticle *particle;
265 RESOLELEC * elementElec;
266 RESOLMUON *elementMuon;
267 RESOLJET *elementJet;
268 TAUHAD *elementTaujet;
269 ETMIS *elementEtmis;
270
271 int numTau=0, numTauRec=0;
272
273 RESOLution *DET = new RESOLution();
274 /*
275 string detectorcard = "data/DetectorCard_CMS.dat";
276 const float dR_jetpairing = 0.25;
277 const float jet_pt_cut = 1;
278 */
279 string detectorcard = "data/DetectorCard_ATLAS.dat";
280 const float dR_jetpairing = 0.2;
281 const float jet_pt_cut = 7;
282 DET->ReadDataCard(detectorcard);
283
284 JetsUtil *JETRUN = new JetsUtil(detectorcard);
285
286 TLorentzVector genMomentum(0,0,0,0); // generator level information
287 TLorentzVector recoMomentum(0,0,0,0);// reconstruction level information
288 LorentzVector jetMomentum;
289
290 vector<fastjet::PseudoJet> input_particlesGEN;//for FastJet algorithm
291 vector<fastjet::PseudoJet> sorted_jetsGEN;
292 vector<int> NTrackJet;
293 vector<TLorentzVector> towers;
294
295 // 2. Loop over all events
296 Long64_t entry, allEntries = treeReader->GetEntries();
297 cout << endl << endl;
298 cout <<"*********************************************************************"<< endl;
299 cout <<"** **"<< endl;
300 cout <<"** **"<< endl;
301 cout <<"** **"<< endl;
302 cout <<"** ####### Start resolution processing ######## **"<< endl;
303 cout << left << setw(52) <<"** Total number of events to run: "<<""
304 << left << setw(15) << allEntries <<""
305 << right << setw(2) <<"**"<<endl;
306
307 ExRootProgressBar *Progress = new ExRootProgressBar(allEntries);
308
309 for(entry = 0; entry < allEntries; ++entry)
310 {
311 Progress->Update(entry);
312 TLorentzVector PTmisReco(0,0,0,0);
313 TLorentzVector PTmisGEN(0,0,0,0);
314 treeReader->ReadEntry(entry);
315 treeReaderGEN->ReadEntry(entry);
316 treeWriter->Clear();
317
318 TSimpleArray<TRootGenParticle> bGen;
319 itGen.Reset();
320 TSimpleArray<TRootGenParticle> NFCentralQ;
321
322 input_particlesGEN.clear();
323 towers.clear();
324
325 // Loop over all particles in event
326 while( (particle = (TRootGenParticle*) itGen.Next()) )
327 {
328 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
329
330 int pid = abs(particle->PID);
331 float eta = fabs(particle->Eta);
332
333 //input generator level particle for jet algorithm
334 if(particle->Status == 1 && eta < DET->CEN_max_calo_fwd)
335 {
336 input_particlesGEN.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
337 }
338
339 //Calculate ETMIS from generated particles
340 if((pid == pNU1) || (pid == pNU2) || (pid == pNU3))PTmisGEN = PTmisGEN + genMomentum;
341
342 //Electrons and muons
343 if( (particle->Status == 1) &&
344 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
345 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
346 )
347 {
348 eta=fabs(genMomentum.Eta());
349
350 switch(pid) {
351
352 case pE: // all electrons with eta < DET->MAX_CALO_FWD
353 PairingElec(recoMomentum,genMomentum,branchElec);
354 if(recoMomentum.E()!=0){
355 elementElec=(RESOLELEC*) branchelec->NewEntry();
356 elementElec->E = genMomentum.E();
357 elementElec->SmearedE = recoMomentum.E();}
358 break; // case pE
359 case pMU: // all muons with eta < DET->MAX_MU
360 PairingMuon(recoMomentum,genMomentum,branchMuon);
361 if(recoMomentum.E()!=0){
362 elementMuon = (RESOLMUON*) branchmuon->NewEntry();
363 elementMuon->OverPT = 1./genMomentum.Pt();
364 elementMuon->OverSmearedPT = 1./recoMomentum.Pt();}
365 break; // case pMU
366 default:
367 break;
368 } // switch (pid)
369 }
370
371 } // while
372
373 //compute missing transverse energy from calo towers
374 TIter itCalo((TCollection*)branchTowers);
375 TRootCalo *calo;
376 itCalo.Reset();
377 TLorentzVector Att(0.,0.,0.,0.);
378 float ScalarEt=0;
379 while( (calo = (TRootCalo*) itCalo.Next()) )
380 {
381 if(calo->E !=0){
382 Att.SetPtEtaPhiE(calo->getET(),calo->Eta,calo->Phi,calo->E);
383 towers.push_back(Att);
384 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
385 {
386 ScalarEt = ScalarEt + calo->getET();
387 PTmisReco = PTmisReco + Att;
388 }
389 }
390 }
391 elementEtmis= (ETMIS*) branchetmis->NewEntry();
392 elementEtmis->Et = (PTmisGEN).Pt();
393 elementEtmis->Ex = (-PTmisGEN).Px();
394 elementEtmis->SEt = ScalarEt;
395 elementEtmis->EtSmeare = (PTmisReco).Pt()-(PTmisGEN).Pt();
396 elementEtmis->ExSmeare = (-PTmisReco).Px()-(PTmisGEN).Px();
397
398 //*****************************
399 sorted_jetsGEN=JETRUN->RunJetsResol(input_particlesGEN);
400
401 TSimpleArray<TRootGenParticle> TausHadr = TauHadr(branchGen);
402 TLorentzVector JETreco(0,0,0,0);
403 for (unsigned int i = 0; i < sorted_jetsGEN.size(); i++) {
404 TLorentzVector JETgen(0,0,0,0);
405 JETgen.SetPxPyPzE(sorted_jetsGEN[i].px(),sorted_jetsGEN[i].py(),sorted_jetsGEN[i].pz(),sorted_jetsGEN[i].E());
406 PairingJet(JETreco,JETgen,branchJet,dR_jetpairing);
407 if(JETreco.Pt()>jet_pt_cut)
408 {
409 elementJet= (RESOLJET*) branchjet->NewEntry();
410 elementJet->PT = JETgen.Et();
411 elementJet->SmearedPT = JETreco.Et()/JETgen.Et();
412 elementJet->E = JETgen.E();
413 elementJet->dE = (JETreco.E()-JETgen.E())/JETgen.E() ;
414 elementJet->dE2 = pow( (JETreco.E()-JETgen.E())/JETgen.E() , 2.) ;
415 }
416 }
417 numTau += TausHadr.GetEntries();
418
419 TIter itJet((TCollection*)branchJet);
420 TRootJet *jet;
421 itJet.Reset();
422 while( (jet = (TRootJet*) itJet.Next()) )
423 {
424 TLorentzVector JETT(0,0,0,0);
425 JETT.SetPxPyPzE(jet->Px,jet->Py,jet->Pz,jet->E);
426 if(fabs(JETT.Eta()) < (DET->CEN_max_tracker - DET->TAU_track_scone))
427 {
428 for(Int_t i=0; i<TausHadr.GetEntries();i++)
429 {
430 if(DeltaR(TausHadr[i]->Phi,TausHadr[i]->Eta,JETT.Phi(),JETT.Eta())<0.1)
431 {
432 elementTaujet= (TAUHAD*) branchtaujet->NewEntry();
433 elementTaujet->EnergieCen = EnergySmallCone(towers,JETT.Eta(),JETT.Phi(),DET->TAU_energy_scone,DET->JET_seed)/JETT.E();
434 elementTaujet->NumTrack = NumTracks(branchTracks,DET->TAU_track_pt,JETT.Eta(),JETT.Phi(),DET->TAU_track_scone);
435 if( (EnergySmallCone(towers,JETT.Eta(),JETT.Phi(),DET->TAU_energy_scone,DET->JET_seed)/JETT.E()) > 0.95
436 && (NumTracks(branchTracks,DET->TAU_track_pt,JETT.Eta(),JETT.Phi(),DET->TAU_track_scone))==1)numTauRec++;
437 }
438 }
439 }
440
441
442 } // for itJet : loop on all jets
443
444 treeWriter->Fill();
445 } // Loop over all events
446 treeWriter->Write();
447 globalwatch.Stop();
448
449
450 // Screen output
451 cout <<"** **"<< endl;
452 cout <<"** **"<< endl;
453 cout <<"** ################## Time report ################# **"<< endl;
454 cout << left << setw(32) <<"** Time report for "<<""
455 << left << setw(15) << allEntries <<""
456 << right << setw(22) <<"events **"<<endl;
457 cout <<"** **"<< endl;
458 cout << left << setw(10) <<"**"<<""
459 << left << setw(15) <<"Time (s):"<<""
460 << right << setw(15) <<"CPU"<<""
461 << right << setw(15) <<"Real"<<""
462 << right << setw(14) <<"**"<<endl;
463 cout << left << setw(10) <<"**"<<""
464 << left << setw(15) <<" + Global:"<<""
465 << right << setw(15) <<globalwatch.CpuTime()<<""
466 << right << setw(15) <<globalwatch.RealTime()<<""
467 << right << setw(14) <<"**"<<endl;
468
469
470 string tempstring = detectorcard + " has been used.";
471 cout <<"** ################## Configuration ############### **"<< endl;
472 cout << left << setw(16) << "** " << ""
473 << left << setw(51) << tempstring << "**" << endl;
474
475 tempstring = outputfilename + " has been created.";
476 cout << left << setw(16) << "** " << ""
477 << left << setw(51) << tempstring << "**" << endl;
478
479 cout << left << setw(16) << "** " << ""
480 << left << setw(51) << get_time_date() << "**" << endl;
481
482 cout <<"** **"<< endl;
483 cout <<"** Exiting ... **"<< endl;
484 cout <<"** **"<< endl;
485 cout <<"** **"<< endl;
486 cout <<"*********************************************************************"<< endl;
487
488
489 delete treeWriter;
490 delete treeReader;
491 delete DET;
492}
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
Note: See TracBrowser for help on using the repository browser.