Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 566

Last change on this file since 566 was 564, checked in by Xavier Rouby, 14 years ago

a muon/electron is not reconstruced if no associated track

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