Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 368

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

update V1.6

File size: 36.7 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \----------------------------------------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
12** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
13** FROG: [hep-ex/0901.2718v1] **
14** **
15** ------------------------------------------------------------------ **
16** **
17** Main authors: **
18** ------------- **
19** **
20** Severine Ovyn Xavier Rouby **
21** severine.ovyn@uclouvain.be xavier.rouby@cern **
22** **
23** Center for Particle Physics and Phenomenology (CP3) **
24** Universite catholique de Louvain (UCL) **
25** Louvain-la-Neuve, Belgium **
26** **
27** Copyright (C) 2008-2009, **
28** All rights reserved. **
29** **
30***********************************************************************/
31
32/// \file Delphes.cpp
33/// \brief Executable for Delphes
34
35#include "TChain.h"
36#include "TApplication.h"
37#include "TStopwatch.h"
38#include "TFile.h"
39
40#include "ExRootTreeReader.h"
41#include "ExRootTreeWriter.h"
42#include "ExRootTreeBranch.h"
43#include "ExRootProgressBar.h"
44
45#include "DataConverter.h"
46#include "LHEFConverter.h"
47#include "HepMCConverter.h"
48#include "HEPEVTConverter.h"
49#include "STDHEPConverter.h"
50#include "LHCOConverter.h"
51
52#include "SmearUtil.h"
53#include "CaloUtil.h"
54#include "BFieldProp.h"
55#include "TriggerUtil.h"
56#include "VeryForward.h"
57#include "JetsUtil.h"
58#include "FrogUtil.h"
59
60#include <vector>
61#include <iostream>
62
63using namespace std;
64
65//------------------------------------------------------------------------------
66void todo(string filename) {
67 ifstream infile(filename.c_str());
68 cout << "** TODO list ..." << endl;
69 while(infile.good()) {
70 string temp;
71 getline(infile,temp);
72 cout << "*" << temp << endl;
73 }
74 cout << "** done...\n";
75}
76
77//------------------------------------------------------------------------------
78
79int main(int argc, char *argv[])
80{
81
82 int appargc = 2;
83 char *appName= new char[20];
84 char *appOpt= new char[20];
85 sprintf(appName,"Delphes");
86 sprintf(appOpt,"-b");
87 char *appargv[] = {appName,appOpt};
88 TApplication app(appName, &appargc, appargv);
89 delete [] appName;
90 delete [] appOpt;
91
92 if(argc != 3 && argc != 4 && argc != 5) {
93 cout << " Usage: " << argv[0] << " input_file output_file [detector_card] [trigger_card] " << endl;
94 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
95 cout << " output_file - output file." << endl;
96 cout << " detector_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
97 cout << " trigger_card - Datacard containing the trigger algorithms (optional) "<<endl;
98 exit(1);
99 }
100
101
102
103
104 cout << endl << endl;
105
106 cout <<"*********************************************************************"<< endl;
107 cout <<"*********************************************************************"<< endl;
108 cout <<"** **"<< endl;
109 cout <<"** Welcome to **"<< endl;
110 cout <<"** **"<< endl;
111 cout <<"** **"<< endl;
112 cout <<"** .ddddddd- lL hH **"<< endl;
113 cout <<"** -Dd` `dD: Ll hH` **"<< endl;
114 cout <<"** dDd dDd eeee. lL .pp+pp Hh+hhh` -eeee- `sssss **"<< endl;
115 cout <<"** -Dd `DD ee. ee Ll .Pp. PP Hh. HH. ee. ee sSs **"<< endl;
116 cout <<"** dD` dDd eEeee: lL. pP. pP hH hH` eEeee:` -sSSSs. **"<< endl;
117 cout <<"** .Dd :dd eE. LlL PpppPP Hh Hh eE sSS **"<< endl;
118 cout <<"** dddddd:. eee+: lL. pp. hh. hh eee+ sssssS **"<< endl;
119 cout <<"** Pp **"<< endl;
120 cout <<"** **"<< endl;
121 cout <<"** Delphes, a framework for the fast simulation **"<< endl;
122 cout <<"** of a generic collider experiment **"<< endl;
123 cout <<"** arXiv:0903.2225v1 [hep-ph] **"<< endl;
124 cout <<"** **"<< endl;
125 cout <<"** --- Version 1.6 of Delphes --- **"<< endl;
126 cout <<"** Last date of change: 7 May 2009 **"<< endl;
127 cout <<"** **"<< endl;
128 cout <<"** **"<< endl;
129 cout <<"** This package uses: **"<< endl;
130 cout <<"** ------------------ **"<< endl;
131 cout <<"** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **"<< endl;
132 cout <<"** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **"<< endl;
133 cout <<"** FROG: [hep-ex/0901.2718v1] **"<< endl;
134 cout <<"** **"<< endl;
135 cout <<"**-----------------------------------------------------------------**"<< endl;
136 cout <<"** **"<< endl;
137 cout <<"** Main authors: **"<< endl;
138 cout <<"** ------------- **"<< endl;
139 cout <<"** **"<< endl;
140 cout <<"** Séverine Ovyn Xavier Rouby **"<< endl;
141 cout <<"** severine.ovyn@uclouvain.be xavier.rouby@cern **"<< endl;
142 cout <<"** Center for Particle Physics and Phenomenology (CP3) **"<< endl;
143 cout <<"** Universite Catholique de Louvain (UCL) **"<< endl;
144 cout <<"** Louvain-la-Neuve, Belgium **"<< endl;
145 cout <<"** **"<< endl;
146 cout <<"**-----------------------------------------------------------------**"<< endl;
147 cout <<"** **"<< endl;
148 cout <<"** Former Delphes versions and documentation can be found on : **"<< endl;
149 cout <<"** http://www.fynu.ucl.ac.be/delphes.html **"<< endl;
150 cout <<"** **"<< endl;
151 cout <<"** **"<< endl;
152 cout <<"** Disclaimer: this program is a beta version of Delphes and **"<< endl;
153 cout <<"** therefore comes without guarantees. Beware of errors and please **"<< endl;
154 cout <<"** give us your feedbacks about potential bugs **"<< endl;
155 cout <<"** **"<< endl;
156 cout <<"*********************************************************************"<< endl;
157 cout <<"*********************************************************************"<< endl;
158
159 // 1. ********** initialisation ***********
160
161 srand (time (NULL)); /* Initialisation du générateur */
162 TStopwatch globalwatch, loopwatch, triggerwatch, frogwatch, lhcowatch;
163 globalwatch.Start();
164
165
166 //read the output TROOT file
167 string inputFileList(argv[1]), outputfilename(argv[2]);
168 if(outputfilename.find(".root") > outputfilename.length()) {
169 cout <<"** ERROR: 'output_file' should be a .root file. Exiting... **"<< endl;
170 exit(1);
171 }
172 //create output log-file name
173 string forLog = outputfilename;
174 string LogName = forLog.erase(forLog.find(".root"));
175 LogName = LogName+"_run.log";
176
177 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
178 outputFile->Close();
179
180 string line;
181 ifstream infile(inputFileList.c_str());
182 infile >> line; // the first line determines the type of input files
183
184 //read the datacard input file
185 string DetDatacard("data/DetectorCard.dat"); //for detector smearing parameters
186 string TrigDatacard("data/TriggerCard.dat"); //for trigger selection
187
188 string lineCard1,lineCard2;
189 bool detecCard=false,trigCard=false;
190 if(argv[3])
191 {
192 ifstream infile1(argv[3]);
193 infile1 >> lineCard1; // the first line determines the type of input files
194 if(strstr(lineCard1.c_str(),"DETECTOR") && detecCard==true)
195 cerr <<"** ERROR: A DETECTOR card has already been loaded **"<< endl;
196 else if(strstr(lineCard1.c_str(),"DETECTOR") && detecCard==false){DetDatacard =argv[3]; detecCard=true;}
197 else if(strstr(lineCard1.c_str(),"TRIGGER") && trigCard==true)
198 cerr <<"** ERROR: A TRIGGER card has already been loaded **"<< endl;
199 else if(strstr(lineCard1.c_str(),"TRIGGER") && trigCard==false){TrigDatacard =argv[3]; trigCard=true;}
200 }
201 if(argv[4])
202 {
203 ifstream infile2(argv[4]);
204 infile2 >> lineCard2; // the first line determines the type of input files
205 if(strstr(lineCard2.c_str(),"DETECTOR") && detecCard==true)
206 cerr <<"** ERROR: A DETECTOR card has already been loaded **"<< endl;
207 else if(strstr(lineCard2.c_str(),"DETECTOR") && detecCard==false){DetDatacard =argv[4]; detecCard=true;}
208 else if(strstr(lineCard2.c_str(),"TRIGGER") && trigCard==true)
209 cerr <<"** ERROR: A TRIGGER card has already been loaded **"<< endl;
210 else if(strstr(lineCard2.c_str(),"TRIGGER") && trigCard==false){TrigDatacard =argv[4]; trigCard=true;}
211 }
212
213 //Smearing information
214 RESOLution *DET = new RESOLution();
215 cout <<"** **"<< endl;
216 cout <<"** ####### Start reading DETECTOR parameters ####### **"<< endl;
217 cout << left << setw(40) <<"** Opening configuration card: "<<""
218 << left << setw(27) << DetDatacard <<""
219 << right << setw(2) <<"**"<<""<<endl;
220 DET->ReadDataCard(DetDatacard);
221 DET->Logfile(LogName);
222 cout << left << setw(40) <<"** Parameters summarised in: "<<""
223 << left << setw(27) << LogName <<""
224 << right << setw(2) <<"**"<<""<<endl;
225 cout <<"** **"<< endl;
226
227 //Trigger information
228 cout <<"** ########### Start reading TRIGGER card ########## **"<< endl;
229 if(trigCard==false)
230 {
231 cout <<"** WARNING: Datacard not found, use default card **" << endl;
232 TrigDatacard="data/TriggerCard.dat";
233 }
234 TriggerTable *TRIGT = new TriggerTable();
235 TRIGT->TriggerCardReader(TrigDatacard.c_str());
236 TRIGT->PrintTriggerTable(LogName);
237 if(DET->FLAG_trigger == 1)
238 {
239 cout << left << setw(40) <<"** Opening configuration card: "<<""
240 << left << setw(27) << TrigDatacard <<""
241 << right << setw(2) <<"**"<<""<<endl;
242 cout <<"** **"<< endl;
243 }
244
245 //Propagation of tracks in the B field
246 TrackPropagation *TRACP = new TrackPropagation(DET);
247
248 //Jet information
249 JetsUtil *JETRUN = new JetsUtil(DET);
250
251 //VFD information
252 VeryForward * VFD = new VeryForward(DET);
253
254 // data converters
255 cout <<"** **"<<endl;
256 cout <<"** ####### Start convertion to TRoot format ######## **"<< endl;
257
258 if(line.length() == 1+line.find_last_of(".hep"))
259 {
260 cout <<"** StdHEP file format detected **"<<endl;
261 cout <<"** This can take several minutes **"<< endl;
262 STDHEPConverter converter(inputFileList,outputfilename);//case ntpl file in input list
263 }
264 else if(line.length() == 1+line.find_last_of(".lhe"))
265 {
266 cout <<"** LHEF file format detected **"<<endl;
267 cout <<"** This can take several minutes **"<< endl;
268 LHEFConverter converter(inputFileList,outputfilename);//case ntpl file in input list
269 }
270 else if(line.length() == 1+line.find_last_of(".root"))
271 {
272 cout <<"** h2root file format detected **"<<endl;
273 cout <<"** This can take several minutes **"<< endl;
274 HEPEVTConverter converter(inputFileList,outputfilename);//case ntpl file in input list
275 }
276 else if(line.length() == 1+line.find_last_of(".hepmc"))
277 {
278 cout <<"** HepMC ASCII file format detected **"<<endl;
279 cout <<"** This can take several minutes **"<< endl;
280 HepMCConverter converter(inputFileList,outputfilename);
281 }
282 else {
283 cerr << left << setw(4) <<"** "<<""
284 << left << setw(63) << line.c_str() <<""
285 << right << setw(2) <<"**"<<endl;
286 cerr <<"** ERROR: File format not identified -- Exiting... **"<< endl;
287 cout <<"** **"<< endl;
288 cout <<"*********************************************************************"<< endl;
289 return -1;};
290 cout <<"** Exiting conversion... **"<< endl;
291
292 TChain chain("GEN");
293 chain.Add(outputfilename.c_str());
294 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
295 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
296
297 TIter itGen((TCollection*)branchGen);
298
299 //Output file : contents of the analysis object data
300 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
301 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
302 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
303 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
304 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
305 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
306 ExRootTreeBranch *branchTrack = treeWriter->NewBranch("Tracks", TRootTracks::Class());
307 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
308 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
309 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
310 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
311 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
312
313 TRootETmis *elementEtmis;
314 TRootElectron *elementElec;
315 TRootMuon *elementMu;
316 TRootPhoton *elementPhoton;
317 TRootTracks * elementTrack;
318 TRootCalo *elementCalo;
319
320 TLorentzVector genMomentum(0,0,0,0); // four-momentum at the vertex
321 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
322 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
323 LorentzVector jetMomentum;
324
325 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
326 vector<fastjet::PseudoJet> sorted_jets;
327 vector<TRootTracks> TrackCentral;
328 vector<PhysicsTower> towers;
329 vector<D_Particle> electron;
330 vector<D_Particle> muon;
331 vector<D_Particle> gamma;
332
333 vector<int> NTrackJet;
334
335 TSimpleArray<TRootC::GenParticle> NFCentralQ;
336
337 D_CaloList list_of_calorimeters;
338 D_CaloElement CentralCalo("centralcalo",
339 -DET->CEN_max_calo_cen, DET->CEN_max_calo_cen,
340 DET->ELG_Ccen, DET->ELG_Ncen, DET->ELG_Scen,
341 DET->HAD_Chcal, DET->HAD_Nhcal, DET->HAD_Shcal);
342 D_CaloElement ForwardCalo("forwardcalo",
343 DET->CEN_max_calo_cen, DET->CEN_max_calo_fwd,
344 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
345 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
346 D_CaloElement BackwardCalo("backwardcalo",
347 -DET->CEN_max_calo_fwd, -DET->CEN_max_calo_cen,
348 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
349 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
350 //D_CaloElement CastorCalo("castor",5.5,6.6,1,0,0,1,0,0);
351 list_of_calorimeters.addElement(CentralCalo);
352 list_of_calorimeters.addElement(ForwardCalo);
353 list_of_calorimeters.addElement(BackwardCalo);
354 //list_of_calorimeters.addElement(CastorCalo);
355 list_of_calorimeters.sortElements();
356
357
358 // 2. ********** Loop over all events ***********
359 Long64_t entry, allEntries = treeReader->GetEntries();
360 cout <<"** **"<<endl;
361 cout <<"** ####### Start fast detector simulation ######## **"<< endl;
362 cout << left << setw(52) <<"** Total number of events to run: "<<""
363 << left << setw(15) << allEntries <<""
364 << right << setw(2) <<"**"<<endl;
365
366 ExRootProgressBar *Progress = new ExRootProgressBar(allEntries);
367
368 loopwatch.Start();
369
370 // loop on all events
371 for(entry = 0; entry < allEntries; ++entry)
372 {
373 Progress->Update(entry);
374 TLorentzVector PTmis(0,0,0,0);
375 treeReader->ReadEntry(entry);
376 treeWriter->Clear();
377
378 electron.clear();
379 muon.clear();
380 gamma.clear();
381 NFCentralQ.Clear();
382
383 TrackCentral.clear();
384 towers.clear();
385 input_particles.clear();
386 NTrackJet.clear();
387
388
389 D_CaloTowerList list_of_active_towers;
390 D_CaloTowerList list_of_towers_with_photon; // to speed up the code: will only look in interesting towers for gamma candidates
391
392
393 // 2.1a Loop over all particles in event, to fill the towers
394 itGen.Reset();
395 TRootC::GenParticle *particleG;
396 while( (particleG = (TRootC::GenParticle*) itGen.Next()) )
397 {
398 TRootGenParticle *particle = new TRootGenParticle(particleG);
399 int pid = abs(particle->PID);
400 particle->Charge=ChargeVal(particle->PID);
401 particle->setFractions(); // init
402
403
404 // 2.1a.1********************* preparation for the b-tagging
405 //// This subarray is needed for the B-jet algorithm
406 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
407 if( (pid <= pB || pid == pGLUON) &&// is it a light quark or a gluon, i.e. is it one of these : u,d,c,s,b,g ?
408 fabs(particle->Eta) < DET->CEN_max_tracker &&
409 particle->Status != 1 &&
410 particle->PT > DET->PT_QUARKS_MIN )
411 {
412 NFCentralQ.Add(particleG);
413 }
414
415 // 2.1a.2********************* visible particles only
416 if( (particle->Status == 1) && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) )
417 {
418 // 2.1a.2.1 Central solenoidal magnetic field
419 TRACP->bfield(particle); // fills in particle->EtaCalo et particle->PhiCalo
420 // 2.1a.2.2 Filling the calorimetric towers -- includes also forward detectors ?
421 // first checks if the charged particles reach the calo!
422 if( DET->FLAG_bfield ||
423 particle->Charge==0 ||
424 (!DET->FLAG_bfield && particle->Charge!=0 && particle->PT > DET->TRACK_ptmin))
425 if(
426 (particle->EtaCalo > list_of_calorimeters.getEtamin() ) &&
427 (particle->EtaCalo < list_of_calorimeters.getEtamax() )
428 )
429 {
430 float iEta=UNDEFINED, iPhi=UNDEFINED;
431 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
432 if (iEta != UNDEFINED && iPhi != UNDEFINED)
433 {
434 D_CaloTower tower(iEta,iPhi); // new tower
435 tower.Set_Eem_Ehad_E_ET(particle->E*particle->getFem() , particle->E*particle->getFhad() );
436 list_of_active_towers.addTower(tower);
437 // this list may contain several times the same calotower, as several particles
438 // may leave some energy in the same calotower
439 // After the loop on particles, identical cells in the list should be merged
440 } // iEta and iPhi must be defined
441 }
442
443 // 2.1a.2.3 charged particles in tracker: energy flow
444 // if bfield not simulated, pt should be high enough to be taken into account
445 // it is supposed here that DET->MAX_calo > DET->CEN_max_tracker > DET->CEN_max_mu > 0
446 if( particle->Charge !=0 &&
447 fabs(particle->EtaCalo)< DET->CEN_max_tracker && // stays in the tracker -> track available
448 ( DET->FLAG_bfield ||
449 (!DET->FLAG_bfield && particle->PT > DET->TRACK_ptmin)
450 )
451 )
452 {
453 // 2.1a.2.3.1 Filling the particle properties + smearing
454 // Hypothesis: the final eta/phi are the ones from the generator, thanks to the track reconstruction
455 // This is the EnergyFlow hypothesis
456 particle->SetEtaPhi(particle->Eta,particle->Phi);
457 float sET=UNDEFINED; // smeared ET, computed from the smeared E -> needed for the tracks
458
459 // 2.1a.2.3.2 Muons
460 if (pid == pMU && fabs(particle->EtaCalo)< DET->CEN_max_mu)
461 {
462 TLorentzVector p;
463 float sPT = gRandom->Gaus(particle->PT, DET->MU_SmearPt*particle->PT );
464 if (sPT > 0 && sPT > DET->PTCUT_muon)
465 {
466 p.SetPtEtaPhiE(sPT,particle->Eta,particle->Phi,sPT*cosh(particle->Eta));
467 muon.push_back(D_Particle(p,pMU,particle->EtaCalo,particle->PhiCalo));
468 }
469 sET = (sPT >0)? sPT : 0;
470 }
471 // 2.1a.2.3.3 Electrons
472 else if (pid == pE)
473 {
474 // Finds in which calorimeter the particle has gone, to know its resolution
475
476 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
477 if(currentCalo.getName() == dummyCalo.getName())
478 {
479 cout << "** Warning: the calo coverage behind the tracker is not complete! **" << endl;
480 }
481
482 // final smeared EM energy // electromagnetic fraction F_em =1 for electrons;
483 float sE = currentCalo.getElectromagneticResolution().Smear(particle->E);
484 if (sE>0)
485 {
486 sET = sE/cosh(particle->Eta);
487 // NB: ET is found via the calorimetry and not via the track curvature
488
489 TLorentzVector p;
490 p.SetPtEtaPhiE(sET,particle->Eta,particle->Phi,sE);
491 if (sET > DET->PTCUT_elec)
492 electron.push_back(D_Particle(p,particle->PID,particle->EtaCalo,particle->PhiCalo));
493 }
494 else { sET=0;} // if negative smeared energy -- needed for the tracks
495 }
496 // 2.1a.2.3.4 Other charged particles : smear them for the tracks!
497 else
498 { //other particles
499 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
500 float sEem = currentCalo.getElectromagneticResolution().Smear(particle->E * particle->getFem());
501 float sEhad = currentCalo.getHadronicResolution().Smear(particle->E * particle->getFhad());
502 float sE = ( (sEem>0)? sEem : 0 ) + ( (sEhad>0)? sEhad : 0 );
503 sET = sE/cosh(particle->EtaCalo);
504 }
505
506 // 2.1a.2.3.5 Tracks
507 if( (rand()%100) < DET->TRACK_eff && sET!=0)
508 {
509 elementTrack = (TRootTracks*) branchTrack->NewEntry();
510 elementTrack->Set(particle->Eta, particle->Phi, particle->EtaCalo, particle->PhiCalo, sET, particle->Charge);
511 TrackCentral.push_back(*elementTrack); // tracks at vertex!
512 // TODO!!! apply a smearing on the position of the origin of the track
513 // TODO!!! elementTracks->SetPositionOut(Xout,Yout,Zout);
514 }
515 } // 2.1a.2.3 : if tracker/energy-flow
516 // 2.1a.2.4 Photons
517 // stays in the tracker -> track available -> gamma ID
518 else if( (pid == pGAMMA) && fabs(particle->EtaCalo)< DET->CEN_max_tracker )
519 {
520 float iEta=UNDEFINED, iPhi=UNDEFINED;
521 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
522 D_CaloTower tower(iEta,iPhi);
523 // stores the list of towers where to apply the photon ID algorithm. Just a trick for a faster search
524 list_of_towers_with_photon.addTower(tower);
525 }
526 // 2.1a.2.5 : very forward detectors
527 else
528 {
529 if (DET->FLAG_RP==1)
530 {
531 // for the moment, only protons are transported
532 // BUT !!! could be a beam of other particles! (heavy ions?)
533 // BUT ALSO !!! if very forward muons, or others!
534 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
535 }
536 // 2.1a.2.6: Zero degree calorimeter
537 if(DET->FLAG_vfd==1)
538 {
539 VFD->ZDC(treeWriter,branchZDC,particle);
540 }
541 }
542
543 } // 2.1a.2 : if visible particle
544 delete particle;
545 }// loop on all particles 2.1a
546
547
548 // 2.1b loop on all (activated) towers
549 // at this stage, list_of_active_towers may contain several times the same tower
550 // first step is to merge identical towers, by matching their (iEta,iPhi)
551
552 list_of_active_towers.sortElements();
553 list_of_active_towers.mergeDuplicates();
554
555 // Calotower smearing
556 list_of_active_towers.smearTowers(list_of_calorimeters);
557
558 for(unsigned int i=0; i<list_of_active_towers.size(); i++)
559 {
560 float iEta = list_of_active_towers[i].getEta();
561 float iPhi = list_of_active_towers[i].getPhi();
562 float e = list_of_active_towers[i].getE();
563 if(iEta != UNDEFINED && iPhi != UNDEFINED && e!=0)
564 {
565 elementCalo = (TRootCalo*) branchCalo->NewEntry();
566 elementCalo->set(list_of_active_towers[i]);
567 // not beautiful : should be improved!
568 TLorentzVector p;
569 p.SetPtEtaPhiE(list_of_active_towers[i].getET(), iEta, iPhi, e );
570 PhysicsTower Tower(LorentzVector(p.Px(),p.Py(),p.Pz(),p.E()));
571 towers.push_back(Tower);
572 }
573 } // loop on towers
574
575 // 2.1c photon ID
576 // list_of_towers_with_photon is the list of towers with photon candidates
577 // already smeared !
578 // sorts the vector and smears duplicates
579 list_of_towers_with_photon.mergeDuplicates();
580 for(unsigned int i=0; i<list_of_towers_with_photon.size(); i++) {
581 float eta = list_of_towers_with_photon[i].getEta();
582 float phi = list_of_towers_with_photon[i].getPhi();
583 D_CaloTower cal(list_of_active_towers.getElement(eta,phi));
584 if(cal.getEta() != UNDEFINED && cal.getPhi() != UNDEFINED && cal.getE() > 0)
585 {
586 TLorentzVector p;
587 p.SetPtEtaPhiE(cal.getET(), eta,phi,cal.getE() );
588 if (cal.getET() > DET->PTCUT_gamma) { gamma.push_back(D_Particle(p,pGAMMA,p.Eta(),p.Phi())); }
589 }
590 } // for -- list of photons
591
592 // 2.2 ********** Output preparation & complex objects ***********
593 // 2.2.1 ********************* sorting collections by decreasing pt
594 DET->SortedVector(electron);
595 float iPhiEl=0,iEtaEl=0,ptisoEl=0;
596 for(unsigned int i=0; i < electron.size(); i++)
597 {
598 elementElec = (TRootElectron*) branchElectron->NewEntry();
599 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
600 elementElec->EtaCalo = electron[i].EtaCalo();
601 elementElec->PhiCalo = electron[i].PhiCalo();
602 elementElec->Charge = sign(electron[i].PID());
603 elementElec->IsolFlag = DET->Isolation(electron[i],TrackCentral,DET->ISOL_PT,DET->ISOL_Cone,ptisoEl);
604 elementElec->IsolPt = ptisoEl;
605 DET->BinEtaPhi(elementElec->PhiCalo,elementElec->EtaCalo,iPhiEl,iEtaEl);
606 D_CaloTower calElec(list_of_active_towers.getElement(iEtaEl,iPhiEl));
607 elementElec->EHoverEE = calElec.getEhad()/calElec.getEem();
608 }
609
610 DET->SortedVector(muon);
611 float iPhiMu=0,iEtaMu=0,ptisoMu=0;
612 for(unsigned int i=0; i < muon.size(); i++)
613 {
614 elementMu = (TRootMuon*) branchMuon->NewEntry();
615 elementMu->Charge = sign(muon[i].PID());
616 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
617 elementMu->EtaCalo = muon[i].EtaCalo();
618 elementMu->PhiCalo = muon[i].PhiCalo();
619 elementMu->IsolFlag = DET->Isolation(muon[i],TrackCentral,DET->ISOL_PT,DET->ISOL_Cone,ptisoMu);
620 elementMu->IsolPt = ptisoMu;
621 DET->BinEtaPhi(elementMu->PhiCalo,elementMu->EtaCalo,iPhiMu,iEtaMu);
622 D_CaloTower calMuon(list_of_active_towers.getElement(iEtaMu,iPhiMu));
623 if( calMuon.getEem() !=0 ) elementMu->EHoverEE = calMuon.getEhad()/calMuon.getEem();
624 else elementMu->EHoverEE = UNDEFINED;
625 elementMu->EtRatio = DET->CaloIsolation(muon[i], list_of_active_towers);
626 }
627
628 DET->SortedVector(gamma);
629 for(unsigned int i=0; i < gamma.size(); i++)
630 {
631 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
632 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
633 D_CaloTower calGamma(list_of_active_towers.getElement(gamma[i].EtaCalo(),gamma[i].PhiCalo()));
634 elementPhoton->EHoverEE = calGamma.getEhad()/calGamma.getEem();
635 }
636
637 // 2.2.2 ************* computes the Missing Transverse Momentum
638 TLorentzVector Att(0.,0.,0.,0.);
639 for(unsigned int i=0; i < towers.size(); i++)
640 {
641 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
642 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
643 {
644 PTmis = PTmis + Att;
645 // create a fastjet::PseudoJet with these components and put it onto
646 // back of the input_particles vector
647 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
648 }
649 }
650 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
651 elementEtmis->ET = (PTmis).Pt();
652 elementEtmis->Phi = (-PTmis).Phi();
653 elementEtmis->Px = (-PTmis).Px();
654 elementEtmis->Py = (-PTmis).Py();
655
656 // 2.2.3 ************* jets, B-tag, tau jets
657 vector<int> NTrackJet; //for number of tracks
658 vector<float> EHADEEM; //for energyHad over energyEm
659 sorted_jets=JETRUN->RunJets(input_particles, TrackCentral,NTrackJet,EHADEEM,list_of_active_towers);
660 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ,NTrackJet,EHADEEM);
661 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral,NTrackJet,EHADEEM);
662
663 treeWriter->Fill();
664 } // 2. Loop over all events ('for' loop)
665
666 cout <<"** Exiting detector simulation... **"<< endl;
667
668
669 treeWriter->Write();
670 delete treeWriter;
671 loopwatch.Stop();
672
673
674
675 // 3. ********** Trigger & Frog ***********
676 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
677 triggerwatch.Start();
678 if(DET->FLAG_trigger == 1)
679 {
680 cout <<"** **"<<endl;
681 cout <<"** ########### Start Trigger selection ########### **"<< endl;
682
683 // input
684 TChain chainT("Analysis");
685 chainT.Add(outputfilename.c_str());
686 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
687
688 // output
689 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
690 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
691 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
692 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
693 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
694 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
695
696 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
697 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
698
699
700 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
701 // loop on all entries
702 for(entryT = 0; entryT < allEntriesT; ++entryT) {
703 treeWriterT->Clear();
704 treeReaderT->ReadEntry(entryT);
705 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
706 treeWriterT->Fill();
707 } // loop on all entries
708 cout <<"** Exiting trigger simulation... **"<< endl;
709
710 treeWriterT->Write();
711 delete treeWriterT;
712 delete treeReaderT;
713 } // trigger
714 triggerwatch.Stop();
715
716
717 // 3.2 ************** FROG display
718 frogwatch.Start();
719 if(DET->FLAG_frog == 1) {
720 cout <<"** **"<<endl;
721 cout <<"** ################## Start FROG ################# **"<< endl;
722
723 FrogDisplay *FROG = new FrogDisplay(DET);
724 FROG->BuildEvents(outputfilename);
725 FROG->BuildGeom();
726 delete FROG;
727 cout <<"** Exiting FROG preparation... **"<< endl;
728 }
729 frogwatch.Stop();
730
731 // 3.3 *************** LHCO output
732 lhcowatch.Start();
733 if(DET->FLAG_lhco == 1){
734 cout <<"** **"<<endl;
735 cout <<"** ############ Start LHCO convertion ############ **"<< endl;
736
737 //LHCOConverter *LHCO = new LHCOConverter(outputfilename,LogNameLHCO);
738 LHCOConverter *LHCO = new LHCOConverter(outputfilename,"");
739 LHCO->CopyRunLogFile();
740 LHCO->ConvertExRootAnalysisToLHCO();
741 delete LHCO;
742 cout <<"** Exiting LHCO converstion... **"<< endl;
743 }
744 lhcowatch.Stop();
745
746
747
748 // 4. ********** End & Exit ***********
749
750 globalwatch.Stop();
751 cout <<"** **"<< endl;
752 cout <<"** ################## Time report ################# **"<< endl;
753 cout << left << setw(32) <<"** Time report for "<<""
754 << left << setw(15) << allEntries <<""
755 << right << setw(22) <<"events **"<<endl;
756 cout <<"** **"<< endl;
757 cout << left << setw(10) <<"**"<<""
758 << left << setw(15) <<"Time (s):"<<""
759 << right << setw(15) <<"CPU"<<""
760 << right << setw(15) <<"Real"<<""
761 << right << setw(14) <<"**"<<endl;
762 cout << left << setw(10) <<"**"<<""
763 << left << setw(15) <<" + Global:"<<""
764 << right << setw(15) <<globalwatch.CpuTime()<<""
765 << right << setw(15) <<globalwatch.RealTime()<<""
766 << right << setw(14) <<"**"<<endl;
767 cout << left << setw(10) <<"**"<<""
768 << left << setw(15) <<" + Events:"<<""
769 << right << setw(15) <<loopwatch.CpuTime()<<""
770 << right << setw(15) <<loopwatch.RealTime()<<""
771 << right << setw(14) <<"**"<<endl;
772 if(DET->FLAG_trigger == 1)
773 {
774 cout << left << setw(10) <<"**"<<""
775 << left << setw(15) <<" + Trigger:"<<""
776 << right << setw(15) <<triggerwatch.CpuTime()<<""
777 << right << setw(15) <<triggerwatch.RealTime()<<""
778 << right << setw(14) <<"**"<<endl;
779 }
780 if(DET->FLAG_frog == 1)
781 {
782 cout << left << setw(10) <<"**"<<""
783 << left << setw(15) <<" + Frog:"<<""
784 << right << setw(15) <<frogwatch.CpuTime()<<""
785 << right << setw(15) <<frogwatch.RealTime()<<""
786 << right << setw(14) <<"**"<<endl;
787 }
788 if(DET->FLAG_lhco == 1)
789 {
790 cout << left << setw(10) <<"**"<<""
791 << left << setw(15) <<" + LHCO:"<<""
792 << right << setw(15) <<lhcowatch.CpuTime()<<""
793 << right << setw(15) <<lhcowatch.RealTime()<<""
794 << right << setw(14) <<"**"<<endl;
795 }
796
797
798 cout <<"** **"<< endl;
799 cout <<"** Exiting Delphes ... **"<< endl;
800 cout <<"** **"<< endl;
801 cout <<"*********************************************************************"<< endl;
802 cout <<"*********************************************************************"<< endl;
803
804
805 delete treeReader;
806 delete DET;
807 delete TRIGT;
808 delete TRACP;
809 delete JETRUN;
810 delete VFD;
811
812 // todo("TODO");
813}
Note: See TracBrowser for help on using the repository browser.