Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 309

Last change on this file since 309 was 307, checked in by severine ovyn, 16 years ago

LHCO advanced

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