Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 397

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

time report moved in SmearUtil+logfile

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