Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 359

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

update converters

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