Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 381

Last change on this file since 381 was 380, checked in by Xavier Rouby, 16 years ago

new PDG table

File size: 37.1 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
216 cout <<"** **"<< endl;
217 cout <<"** ####### Start reading DETECTOR parameters ####### **"<< endl;
218 cout << left << setw(40) <<"** Opening configuration card: "<<""
219 << left << setw(27) << DetDatacard <<""
220 << right << setw(2) <<"**"<<""<<endl;
221 DET->ReadDataCard(DetDatacard);
222 DET->Logfile(LogName);
223 cout << left << setw(40) <<"** Parameters summarised in: "<<""
224 << left << setw(27) << LogName <<""
225 << right << setw(2) <<"**"<<""<<endl;
226 cout <<"** **"<< endl;
227 DET->ReadParticleDataGroupTable();
228 //DET->PDGtable.print();
229
230 //Trigger information
231 cout <<"** ########### Start reading TRIGGER card ########## **"<< endl;
232 if(trigCard==false)
233 {
234 cout <<"** WARNING: Datacard not found, use default card **" << endl;
235 TrigDatacard="data/TriggerCard.dat";
236 }
237 TriggerTable *TRIGT = new TriggerTable();
238 TRIGT->TriggerCardReader(TrigDatacard.c_str());
239 TRIGT->PrintTriggerTable(LogName);
240 if(DET->FLAG_trigger == 1)
241 {
242 cout << left << setw(40) <<"** Opening configuration card: "<<""
243 << left << setw(27) << TrigDatacard <<""
244 << right << setw(2) <<"**"<<""<<endl;
245 cout <<"** **"<< endl;
246 }
247
248 //Propagation of tracks in the B field
249 TrackPropagation *TRACP = new TrackPropagation(DET);
250
251 //Jet information
252 JetsUtil *JETRUN = new JetsUtil(DET);
253
254 //VFD information
255 VeryForward * VFD = new VeryForward(DET);
256
257 // data converters
258 cout <<"** **"<<endl;
259 cout <<"** ####### Start convertion to TRoot format ######## **"<< endl;
260
261 if(line.length() == 1+line.find_last_of(".hep"))
262 {
263 cout <<"** StdHEP file format detected **"<<endl;
264 cout <<"** This can take several minutes **"<< endl;
265 STDHEPConverter converter(inputFileList,outputfilename,DET->PDGtable);//case ntpl file in input list
266 }
267 else if(line.length() == 1+line.find_last_of(".lhe"))
268 {
269 cout <<"** LHEF file format detected **"<<endl;
270 cout <<"** This can take several minutes **"<< endl;
271 LHEFConverter converter(inputFileList,outputfilename,DET->PDGtable);//case ntpl file in input list
272 }
273 else if(line.length() == 1+line.find_last_of(".root"))
274 {
275 cout <<"** h2root file format detected **"<<endl;
276 cout <<"** This can take several minutes **"<< endl;
277 HEPEVTConverter converter(inputFileList,outputfilename,DET->PDGtable);//case ntpl file in input list
278 }
279 else if(line.length() == 1+line.find_last_of(".hepmc"))
280 {
281 cout <<"** HepMC ASCII file format detected **"<<endl;
282 cout <<"** This can take several minutes **"<< endl;
283 HepMCConverter converter(inputFileList,outputfilename,DET->PDGtable);
284 }
285 else {
286 cerr << left << setw(4) <<"** "<<""
287 << left << setw(63) << line.c_str() <<""
288 << right << setw(2) <<"**"<<endl;
289 cerr <<"** ERROR: File format not identified -- Exiting... **"<< endl;
290 cout <<"** **"<< endl;
291 cout <<"*********************************************************************"<< endl;
292 return -1;};
293 cout <<"** Exiting conversion... **"<< endl;
294
295 TChain chain("GEN");
296 chain.Add(outputfilename.c_str());
297 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
298 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
299
300 TIter itGen((TCollection*)branchGen);
301
302 //Output file : contents of the analysis object data
303 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
304 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
305 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
306 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
307 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
308 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
309 ExRootTreeBranch *branchTrack = treeWriter->NewBranch("Tracks", TRootTracks::Class());
310 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
311 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
312 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
313 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
314 //ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootForwardTaggerHits::Class());
315 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
316
317 TRootETmis *elementEtmis;
318 TRootElectron *elementElec;
319 TRootMuon *elementMu;
320 TRootPhoton *elementPhoton;
321 TRootTracks * elementTrack;
322 TRootCalo *elementCalo;
323
324 TLorentzVector genMomentum(0,0,0,0); // four-momentum at the vertex
325 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
326 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
327 LorentzVector jetMomentum;
328
329 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
330 vector<fastjet::PseudoJet> sorted_jets;
331 vector<TRootTracks> TrackCentral;
332 vector<PhysicsTower> towers;
333 vector<D_Particle> electron;
334 vector<D_Particle> muon;
335 vector<D_Particle> gamma;
336
337 vector<int> NTrackJet;
338
339 TSimpleArray<TRootC::GenParticle> NFCentralQ;
340
341 D_CaloList list_of_calorimeters;
342 D_CaloElement CentralCalo("centralcalo",
343 -DET->CEN_max_calo_cen, DET->CEN_max_calo_cen,
344 DET->ELG_Ccen, DET->ELG_Ncen, DET->ELG_Scen,
345 DET->HAD_Chcal, DET->HAD_Nhcal, DET->HAD_Shcal);
346 D_CaloElement ForwardCalo("forwardcalo",
347 DET->CEN_max_calo_cen, DET->CEN_max_calo_fwd,
348 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
349 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
350 D_CaloElement BackwardCalo("backwardcalo",
351 -DET->CEN_max_calo_fwd, -DET->CEN_max_calo_cen,
352 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
353 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
354 //D_CaloElement CastorCalo("castor",5.5,6.6,1,0,0,1,0,0);
355 list_of_calorimeters.addElement(CentralCalo);
356 list_of_calorimeters.addElement(ForwardCalo);
357 list_of_calorimeters.addElement(BackwardCalo);
358 //list_of_calorimeters.addElement(CastorCalo);
359 list_of_calorimeters.sortElements();
360
361
362 // 2. ********** Loop over all events ***********
363 Long64_t entry, allEntries = treeReader->GetEntries();
364 cout <<"** **"<<endl;
365 cout <<"** ####### Start fast detector simulation ######## **"<< endl;
366 cout << left << setw(52) <<"** Total number of events to run: "<<""
367 << left << setw(15) << allEntries <<""
368 << right << setw(2) <<"**"<<endl;
369
370 ExRootProgressBar *Progress = new ExRootProgressBar(allEntries);
371
372 loopwatch.Start();
373
374 // loop on all events
375 for(entry = 0; entry < allEntries; ++entry)
376 {
377 Progress->Update(entry);
378 TLorentzVector PTmis(0,0,0,0);
379 treeReader->ReadEntry(entry);
380 treeWriter->Clear();
381
382 electron.clear();
383 muon.clear();
384 gamma.clear();
385 NFCentralQ.Clear();
386
387 TrackCentral.clear();
388 towers.clear();
389 input_particles.clear();
390 NTrackJet.clear();
391
392
393 D_CaloTowerList list_of_active_towers;
394 D_CaloTowerList list_of_towers_with_photon; // to speed up the code: will only look in interesting towers for gamma candidates
395
396
397 // 2.1a Loop over all particles in event, to fill the towers
398 itGen.Reset();
399 TRootC::GenParticle *particleG;
400 while( (particleG = (TRootC::GenParticle*) itGen.Next()) )
401 {
402 TRootGenParticle *particle = new TRootGenParticle(particleG);
403 PdgParticle pdg_part = DET->PDGtable[particle->PID];
404 particle->Charge = pdg_part.charge();
405 particle->M = pdg_part.mass();
406 //particle->Charge=ChargeVal(particle->PID);
407 particle->setFractions(); // init
408 int pid = abs(particle->PID);
409
410
411 // 2.1a.1********************* preparation for the b-tagging
412 //// This subarray is needed for the B-jet algorithm
413 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
414 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 ?
415 fabs(particle->Eta) < DET->CEN_max_tracker &&
416 particle->Status != 1 &&
417 particle->PT > DET->PT_QUARKS_MIN )
418 {
419 NFCentralQ.Add(particleG);
420 }
421
422 // 2.1a.2********************* visible particles only
423 if( (particle->Status == 1) && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) )
424 {
425 // 2.1a.2.1 Central solenoidal magnetic field
426 TRACP->bfield(particle); // fills in particle->EtaCalo et particle->PhiCalo
427 // 2.1a.2.2 Filling the calorimetric towers -- includes also forward detectors ?
428 // first checks if the charged particles reach the calo!
429 if( DET->FLAG_bfield ||
430 particle->Charge==0 ||
431 (!DET->FLAG_bfield && particle->Charge!=0 && particle->PT > DET->TRACK_ptmin))
432 if(
433 (particle->EtaCalo > list_of_calorimeters.getEtamin() ) &&
434 (particle->EtaCalo < list_of_calorimeters.getEtamax() )
435 )
436 {
437 float iEta=UNDEFINED, iPhi=UNDEFINED;
438 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
439 if (iEta != UNDEFINED && iPhi != UNDEFINED)
440 {
441 D_CaloTower tower(iEta,iPhi); // new tower
442 tower.Set_Eem_Ehad_E_ET(particle->E*particle->getFem() , particle->E*particle->getFhad() );
443 list_of_active_towers.addTower(tower);
444 // this list may contain several times the same calotower, as several particles
445 // may leave some energy in the same calotower
446 // After the loop on particles, identical cells in the list should be merged
447 } // iEta and iPhi must be defined
448 }
449
450 // 2.1a.2.3 charged particles in tracker: energy flow
451 // if bfield not simulated, pt should be high enough to be taken into account
452 // it is supposed here that DET->MAX_calo > DET->CEN_max_tracker > DET->CEN_max_mu > 0
453 if( particle->Charge !=0 &&
454 fabs(particle->EtaCalo)< DET->CEN_max_tracker && // stays in the tracker -> track available
455 ( DET->FLAG_bfield ||
456 (!DET->FLAG_bfield && particle->PT > DET->TRACK_ptmin)
457 )
458 )
459 {
460 // 2.1a.2.3.1 Filling the particle properties + smearing
461 // Hypothesis: the final eta/phi are the ones from the generator, thanks to the track reconstruction
462 // This is the EnergyFlow hypothesis
463 particle->SetEtaPhi(particle->Eta,particle->Phi);
464 float sET=UNDEFINED; // smeared ET, computed from the smeared E -> needed for the tracks
465
466 // 2.1a.2.3.2 Muons
467 if (pid == pMU && fabs(particle->EtaCalo)< DET->CEN_max_mu)
468 {
469 TLorentzVector p;
470 float sPT = gRandom->Gaus(particle->PT, DET->MU_SmearPt*particle->PT );
471 if (sPT > 0 && sPT > DET->PTCUT_muon)
472 {
473 p.SetPtEtaPhiE(sPT,particle->Eta,particle->Phi,sPT*cosh(particle->Eta));
474 muon.push_back(D_Particle(p,pMU,particle->EtaCalo,particle->PhiCalo));
475 }
476 sET = (sPT >0)? sPT : 0;
477 }
478 // 2.1a.2.3.3 Electrons
479 else if (pid == pE)
480 {
481 // Finds in which calorimeter the particle has gone, to know its resolution
482
483 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
484 if(currentCalo.getName() == dummyCalo.getName())
485 {
486 cout << "** Warning: the calo coverage behind the tracker is not complete! **" << endl;
487 }
488
489 // final smeared EM energy // electromagnetic fraction F_em =1 for electrons;
490 float sE = currentCalo.getElectromagneticResolution().Smear(particle->E);
491 if (sE>0)
492 {
493 sET = sE/cosh(particle->Eta);
494 // NB: ET is found via the calorimetry and not via the track curvature
495
496 TLorentzVector p;
497 p.SetPtEtaPhiE(sET,particle->Eta,particle->Phi,sE);
498 if (sET > DET->PTCUT_elec)
499 electron.push_back(D_Particle(p,particle->PID,particle->EtaCalo,particle->PhiCalo));
500 }
501 else { sET=0;} // if negative smeared energy -- needed for the tracks
502 }
503 // 2.1a.2.3.4 Other charged particles : smear them for the tracks!
504 else
505 { //other particles
506 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
507 float sEem = currentCalo.getElectromagneticResolution().Smear(particle->E * particle->getFem());
508 float sEhad = currentCalo.getHadronicResolution().Smear(particle->E * particle->getFhad());
509 float sE = ( (sEem>0)? sEem : 0 ) + ( (sEhad>0)? sEhad : 0 );
510 sET = sE/cosh(particle->EtaCalo);
511 }
512
513 // 2.1a.2.3.5 Tracks
514 if( (rand()%100) < DET->TRACK_eff && sET!=0)
515 {
516 elementTrack = (TRootTracks*) branchTrack->NewEntry();
517 elementTrack->Set(particle->Eta, particle->Phi, particle->EtaCalo, particle->PhiCalo, sET, particle->Charge);
518 TrackCentral.push_back(*elementTrack); // tracks at vertex!
519 // TODO!!! apply a smearing on the position of the origin of the track
520 // TODO!!! elementTracks->SetPositionOut(Xout,Yout,Zout);
521 }
522 } // 2.1a.2.3 : if tracker/energy-flow
523 // 2.1a.2.4 Photons
524 // stays in the tracker -> track available -> gamma ID
525 else if( (pid == pGAMMA) && fabs(particle->EtaCalo)< DET->CEN_max_tracker )
526 {
527 float iEta=UNDEFINED, iPhi=UNDEFINED;
528 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
529 D_CaloTower tower(iEta,iPhi);
530 // stores the list of towers where to apply the photon ID algorithm. Just a trick for a faster search
531 list_of_towers_with_photon.addTower(tower);
532 }
533 // 2.1a.2.5 : very forward detectors
534 else
535 {
536 if (DET->FLAG_RP==1)
537 {
538 // for the moment, only protons are transported
539 // BUT !!! could be a beam of other particles! (heavy ions?)
540 // BUT ALSO !!! if very forward muons, or others!
541 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
542 }
543 // 2.1a.2.6: Zero degree calorimeter
544 if(DET->FLAG_vfd==1)
545 {
546 VFD->ZDC(treeWriter,branchZDC,particle);
547 }
548 }
549
550 } // 2.1a.2 : if visible particle
551 delete particle;
552 }// loop on all particles 2.1a
553
554
555 // 2.1b loop on all (activated) towers
556 // at this stage, list_of_active_towers may contain several times the same tower
557 // first step is to merge identical towers, by matching their (iEta,iPhi)
558
559 list_of_active_towers.sortElements();
560 list_of_active_towers.mergeDuplicates();
561
562 // Calotower smearing
563 list_of_active_towers.smearTowers(list_of_calorimeters);
564
565 for(unsigned int i=0; i<list_of_active_towers.size(); i++)
566 {
567 float iEta = list_of_active_towers[i].getEta();
568 float iPhi = list_of_active_towers[i].getPhi();
569 float e = list_of_active_towers[i].getE();
570 if(iEta != UNDEFINED && iPhi != UNDEFINED && e!=0)
571 {
572 elementCalo = (TRootCalo*) branchCalo->NewEntry();
573 elementCalo->set(list_of_active_towers[i]);
574 // not beautiful : should be improved!
575 TLorentzVector p;
576 p.SetPtEtaPhiE(list_of_active_towers[i].getET(), iEta, iPhi, e );
577 PhysicsTower Tower(LorentzVector(p.Px(),p.Py(),p.Pz(),p.E()));
578 towers.push_back(Tower);
579 }
580 } // loop on towers
581
582 // 2.1c photon ID
583 // list_of_towers_with_photon is the list of towers with photon candidates
584 // already smeared !
585 // sorts the vector and smears duplicates
586 list_of_towers_with_photon.mergeDuplicates();
587 for(unsigned int i=0; i<list_of_towers_with_photon.size(); i++) {
588 float eta = list_of_towers_with_photon[i].getEta();
589 float phi = list_of_towers_with_photon[i].getPhi();
590 D_CaloTower cal(list_of_active_towers.getElement(eta,phi));
591 if(cal.getEta() != UNDEFINED && cal.getPhi() != UNDEFINED && cal.getE() > 0)
592 {
593 TLorentzVector p;
594 p.SetPtEtaPhiE(cal.getET(), eta,phi,cal.getE() );
595 if (cal.getET() > DET->PTCUT_gamma) { gamma.push_back(D_Particle(p,pGAMMA,p.Eta(),p.Phi())); }
596 }
597 } // for -- list of photons
598
599 // 2.2 ********** Output preparation & complex objects ***********
600 // 2.2.1 ********************* sorting collections by decreasing pt
601 DET->SortedVector(electron);
602 float iPhiEl=0,iEtaEl=0,ptisoEl=0;
603 for(unsigned int i=0; i < electron.size(); i++)
604 {
605 elementElec = (TRootElectron*) branchElectron->NewEntry();
606 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
607 elementElec->EtaCalo = electron[i].EtaCalo();
608 elementElec->PhiCalo = electron[i].PhiCalo();
609 elementElec->Charge = sign(electron[i].PID());
610 elementElec->IsolFlag = DET->Isolation(electron[i],TrackCentral,DET->ISOL_PT,DET->ISOL_Cone,ptisoEl);
611 elementElec->IsolPt = ptisoEl;
612 DET->BinEtaPhi(elementElec->PhiCalo,elementElec->EtaCalo,iPhiEl,iEtaEl);
613 D_CaloTower calElec(list_of_active_towers.getElement(iEtaEl,iPhiEl));
614 elementElec->EHoverEE = calElec.getEhad()/calElec.getEem();
615 }
616
617 DET->SortedVector(muon);
618 float iPhiMu=0,iEtaMu=0,ptisoMu=0;
619 for(unsigned int i=0; i < muon.size(); i++)
620 {
621 elementMu = (TRootMuon*) branchMuon->NewEntry();
622 elementMu->Charge = sign(muon[i].PID());
623 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
624 elementMu->EtaCalo = muon[i].EtaCalo();
625 elementMu->PhiCalo = muon[i].PhiCalo();
626 elementMu->IsolFlag = DET->Isolation(muon[i],TrackCentral,DET->ISOL_PT,DET->ISOL_Cone,ptisoMu);
627 elementMu->IsolPt = ptisoMu;
628 DET->BinEtaPhi(elementMu->PhiCalo,elementMu->EtaCalo,iPhiMu,iEtaMu);
629 D_CaloTower calMuon(list_of_active_towers.getElement(iEtaMu,iPhiMu));
630 if( calMuon.getEem() !=0 ) elementMu->EHoverEE = calMuon.getEhad()/calMuon.getEem();
631 else elementMu->EHoverEE = UNDEFINED;
632 elementMu->EtRatio = DET->CaloIsolation(muon[i], list_of_active_towers);
633 }
634
635 DET->SortedVector(gamma);
636 for(unsigned int i=0; i < gamma.size(); i++)
637 {
638 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
639 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
640 D_CaloTower calGamma(list_of_active_towers.getElement(gamma[i].EtaCalo(),gamma[i].PhiCalo()));
641 elementPhoton->EHoverEE = calGamma.getEhad()/calGamma.getEem();
642 }
643
644 // 2.2.2 ************* computes the Missing Transverse Momentum
645 TLorentzVector Att(0.,0.,0.,0.);
646 for(unsigned int i=0; i < towers.size(); i++)
647 {
648 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
649 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
650 {
651 PTmis = PTmis + Att;
652 // create a fastjet::PseudoJet with these components and put it onto
653 // back of the input_particles vector
654 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
655 }
656 }
657 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
658 elementEtmis->ET = (PTmis).Pt();
659 elementEtmis->Phi = (-PTmis).Phi();
660 elementEtmis->Px = (-PTmis).Px();
661 elementEtmis->Py = (-PTmis).Py();
662
663 // 2.2.3 ************* jets, B-tag, tau jets
664 vector<int> NTrackJet; //for number of tracks
665 vector<float> EHADEEM; //for energyHad over energyEm
666 sorted_jets=JETRUN->RunJets(input_particles, TrackCentral,NTrackJet,EHADEEM,list_of_active_towers);
667 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ,NTrackJet,EHADEEM);
668 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral,NTrackJet,EHADEEM);
669
670 treeWriter->Fill();
671 } // 2. Loop over all events ('for' loop)
672
673 cout <<"** Exiting detector simulation... **"<< endl;
674
675
676 treeWriter->Write();
677 delete treeWriter;
678 loopwatch.Stop();
679
680
681
682 // 3. ********** Trigger & Frog ***********
683 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
684 triggerwatch.Start();
685 if(DET->FLAG_trigger == 1)
686 {
687 cout <<"** **"<<endl;
688 cout <<"** ########### Start Trigger selection ########### **"<< endl;
689
690 // input
691 TChain chainT("Analysis");
692 chainT.Add(outputfilename.c_str());
693 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
694
695 // output
696 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
697 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
698 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
699 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
700 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
701 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
702
703 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
704 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
705
706
707 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
708 // loop on all entries
709 for(entryT = 0; entryT < allEntriesT; ++entryT) {
710 treeWriterT->Clear();
711 treeReaderT->ReadEntry(entryT);
712 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
713 treeWriterT->Fill();
714 } // loop on all entries
715 cout <<"** Exiting trigger simulation... **"<< endl;
716
717 treeWriterT->Write();
718 delete treeWriterT;
719 delete treeReaderT;
720 } // trigger
721 triggerwatch.Stop();
722
723
724 // 3.2 ************** FROG display
725 frogwatch.Start();
726 if(DET->FLAG_frog == 1) {
727 cout <<"** **"<<endl;
728 cout <<"** ################## Start FROG ################# **"<< endl;
729
730 FrogDisplay *FROG = new FrogDisplay(DET);
731 FROG->BuildEvents(outputfilename);
732 FROG->BuildGeom();
733 delete FROG;
734 cout <<"** Exiting FROG preparation... **"<< endl;
735 }
736 frogwatch.Stop();
737
738 // 3.3 *************** LHCO output
739 lhcowatch.Start();
740 if(DET->FLAG_lhco == 1){
741 cout <<"** **"<<endl;
742 cout <<"** ############ Start LHCO convertion ############ **"<< endl;
743
744 //LHCOConverter *LHCO = new LHCOConverter(outputfilename,LogNameLHCO);
745 LHCOConverter *LHCO = new LHCOConverter(outputfilename,"");
746 LHCO->CopyRunLogFile();
747 LHCO->ConvertExRootAnalysisToLHCO();
748 delete LHCO;
749 cout <<"** Exiting LHCO converstion... **"<< endl;
750 }
751 lhcowatch.Stop();
752
753
754
755 // 4. ********** End & Exit ***********
756
757 globalwatch.Stop();
758 cout <<"** **"<< endl;
759 cout <<"** ################## Time report ################# **"<< endl;
760 cout << left << setw(32) <<"** Time report for "<<""
761 << left << setw(15) << allEntries <<""
762 << right << setw(22) <<"events **"<<endl;
763 cout <<"** **"<< endl;
764 cout << left << setw(10) <<"**"<<""
765 << left << setw(15) <<"Time (s):"<<""
766 << right << setw(15) <<"CPU"<<""
767 << right << setw(15) <<"Real"<<""
768 << right << setw(14) <<"**"<<endl;
769 cout << left << setw(10) <<"**"<<""
770 << left << setw(15) <<" + Global:"<<""
771 << right << setw(15) <<globalwatch.CpuTime()<<""
772 << right << setw(15) <<globalwatch.RealTime()<<""
773 << right << setw(14) <<"**"<<endl;
774 cout << left << setw(10) <<"**"<<""
775 << left << setw(15) <<" + Events:"<<""
776 << right << setw(15) <<loopwatch.CpuTime()<<""
777 << right << setw(15) <<loopwatch.RealTime()<<""
778 << right << setw(14) <<"**"<<endl;
779 if(DET->FLAG_trigger == 1)
780 {
781 cout << left << setw(10) <<"**"<<""
782 << left << setw(15) <<" + Trigger:"<<""
783 << right << setw(15) <<triggerwatch.CpuTime()<<""
784 << right << setw(15) <<triggerwatch.RealTime()<<""
785 << right << setw(14) <<"**"<<endl;
786 }
787 if(DET->FLAG_frog == 1)
788 {
789 cout << left << setw(10) <<"**"<<""
790 << left << setw(15) <<" + Frog:"<<""
791 << right << setw(15) <<frogwatch.CpuTime()<<""
792 << right << setw(15) <<frogwatch.RealTime()<<""
793 << right << setw(14) <<"**"<<endl;
794 }
795 if(DET->FLAG_lhco == 1)
796 {
797 cout << left << setw(10) <<"**"<<""
798 << left << setw(15) <<" + LHCO:"<<""
799 << right << setw(15) <<lhcowatch.CpuTime()<<""
800 << right << setw(15) <<lhcowatch.RealTime()<<""
801 << right << setw(14) <<"**"<<endl;
802 }
803
804
805 cout <<"** **"<< endl;
806 cout <<"** Exiting Delphes ... **"<< endl;
807 cout <<"** **"<< endl;
808 cout <<"*********************************************************************"<< endl;
809 cout <<"*********************************************************************"<< endl;
810
811
812 delete treeReader;
813 delete DET;
814 delete TRIGT;
815 delete TRACP;
816 delete JETRUN;
817 delete VFD;
818
819 // todo("TODO");
820}
Note: See TracBrowser for help on using the repository browser.