Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 473

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

property 'invisible()' added to PdgParticle

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