Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 425

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

bug removed on bad inputfile identification ; Nevent implemented in DetectorCard

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