Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 269

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

nettoyage

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