Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 293

Last change on this file since 293 was 290, checked in by Xavier Rouby, 16 years ago

tracks and Taujets have charge

File size: 34.8 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<TRootTracks> TrackCentral;
322 vector<PhysicsTower> towers;
323 //vector<ParticleUtil> electron;
324 vector<D_Particle> electron;
325 //vector<ParticleUtil> muon;
326 vector<D_Particle> muon;
327 //vector<ParticleUtil> gamma;
328 vector<D_Particle> gamma;
329
330bool FLAG_lhco = true;
331
332 TSimpleArray<TRootGenParticle> NFCentralQ;
333
334D_CaloList list_of_calorimeters;
335D_CaloElement CentralCalo("centralcalo",
336 -DET->CEN_max_calo_cen, DET->CEN_max_calo_cen,
337 DET->ELG_Ccen, DET->ELG_Ncen, DET->ELG_Scen,
338 DET->HAD_Chcal, DET->HAD_Nhcal, DET->HAD_Shcal);
339D_CaloElement ForwardCalo("forwardcalo",
340 DET->CEN_max_calo_cen, DET->CEN_max_calo_fwd,
341 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
342 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
343D_CaloElement BackwardCalo("backwardcalo",
344 -DET->CEN_max_calo_fwd, -DET->CEN_max_calo_cen,
345 DET->ELG_Cfwd, DET->ELG_Nfwd, DET->ELG_Sfwd,
346 DET->HAD_Chf, DET->HAD_Nhf, DET->HAD_Shf );
347//D_CaloElement CastorCalo("castor",5.5,6.6,1,0,0,1,0,0);
348list_of_calorimeters.addElement(CentralCalo);
349list_of_calorimeters.addElement(ForwardCalo);
350list_of_calorimeters.addElement(BackwardCalo);
351//list_of_calorimeters.addElement(CastorCalo);
352list_of_calorimeters.sortElements();
353
354
355// 2. ********** Loop over all events ***********
356 Long64_t entry, allEntries = treeReader->GetEntries();
357 cout <<"** **"<<endl;
358 cout <<"** ####### Start fast detector simulation ######## **"<< endl;
359 cout << left << setw(52) <<"** Total number of events to run: "<<""
360 << left << setw(15) << allEntries <<""
361 << right << setw(2) <<"**"<<endl;
362
363 ExRootProgressBar *Progress = new ExRootProgressBar(allEntries);
364
365 loopwatch.Start();
366
367 // loop on all events
368 for(entry = 0; entry < allEntries; ++entry)
369 {
370 Progress->Update(entry);
371 TLorentzVector PTmis(0,0,0,0);
372 treeReader->ReadEntry(entry);
373 treeWriter->Clear();
374
375 electron.clear();
376 muon.clear();
377 gamma.clear();
378 NFCentralQ.Clear();
379
380 TrackCentral.clear();
381 towers.clear();
382 input_particles.clear();
383
384 {
385 D_CaloTowerList list_of_active_towers;
386 D_CaloTowerList list_of_towers_with_photon; // to speed up the code: will only look in interesting towers for gamma candidates
387
388 // 2.1a Loop over all particles in event, to fill the towers
389 itGen.Reset();
390 GenParticle *particleG;
391 while( (particleG = (GenParticle*) itGen.Next()) ) {
392
393 TRootGenParticle *particle = new TRootGenParticle(particleG);
394 int pid = abs(particle->PID);
395 particle->Charge=ChargeVal(particle->PID);
396 particle->setFractions(); // init
397
398
399 // 2.1a.1********************* preparation for the b-tagging
400 //// This subarray is needed for the B-jet algorithm
401 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
402 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 ?
403 fabs(particle->Eta) < DET->CEN_max_tracker &&
404 particle->Status != 1 &&
405 particle->PT > DET->PT_QUARKS_MIN ) {
406 NFCentralQ.Add(particle);
407 }
408
409 // 2.1a.2********************* visible particles only
410 if( (particle->Status == 1) && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) ){
411
412 // 2.1a.2.1 Central solenoidal magnetic field
413 TRACP->bfield(particle); // fills in particle->EtaCalo et particle->PhiCalo
414 // 2.1a.2.2 Filling the calorimetric towers -- includes also forward detectors ?
415 // first checks if the charged particles reach the calo!
416 if( DET->FLAG_bfield ||
417 particle->Charge==0 ||
418 (!DET->FLAG_bfield && particle->Charge!=0 && particle->PT > DET->TRACK_ptmin))
419 if(
420 (particle->EtaCalo > list_of_calorimeters.getEtamin() ) &&
421 (particle->EtaCalo < list_of_calorimeters.getEtamax() )
422 ) {
423 float iEta=UNDEFINED, iPhi=UNDEFINED;
424 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
425 if (iEta != UNDEFINED && iPhi != UNDEFINED)
426 {
427 D_CaloTower tower(iEta,iPhi); // new tower
428 tower.Set_Eem_Ehad_E_ET(particle->E*particle->getFem() , particle->E*particle->getFhad() );
429
430 list_of_active_towers.addTower(tower);
431 // this list may contain several times the same calotower, as several particles
432 // may leave some energy in the same calotower
433 // After the loop on particles, identical cells in the list should be merged
434 } // iEta and iPhi must be defined
435 }
436
437 // 2.1a.2.3 charged particles in tracker: energy flow
438 // if bfield not simulated, pt should be high enough to be taken into account
439 // it is supposed here that DET->MAX_calo > DET->CEN_max_tracker > DET->CEN_max_mu > 0
440 if( particle->Charge !=0 &&
441 fabs(particle->EtaCalo)< DET->CEN_max_tracker && // stays in the tracker -> track available
442 ( DET->FLAG_bfield ||
443 (!DET->FLAG_bfield && particle->PT > DET->TRACK_ptmin)
444 )
445 ) {
446 // 2.1a.2.3.1 Filling the particle properties + smearing
447 // Hypothesis: the final eta/phi are the ones from the generator, thanks to the track reconstruction
448 // This is the EnergyFlow hypothesis
449 particle->SetEtaPhi(particle->Eta,particle->Phi);
450 float sET=UNDEFINED; // smeared ET, computed from the smeared E -> needed for the tracks
451
452 // 2.1a.2.3.2 Muons
453 if (pid == pMU && fabs(particle->EtaCalo)< DET->CEN_max_mu) {
454 TLorentzVector p;
455 float sPT = gRandom->Gaus(particle->PT, DET->MU_SmearPt*particle->PT );
456 if (sPT > 0 && sPT > DET->PTCUT_muon) {
457 p.SetPtEtaPhiE(sPT,particle->Eta,particle->Phi,sPT*cosh(particle->Eta));
458 muon.push_back(D_Particle(p,pMU,particle->EtaCalo,particle->PhiCalo));
459 }
460 sET = (sPT >0)? sPT : 0;
461 }
462 // 2.1a.2.3.3 Electrons
463 else if (pid == pE) {
464 // Finds in which calorimeter the particle has gone, to know its resolution
465
466 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
467 if(currentCalo.getName() == dummyCalo.getName()) {
468 cout << "** Warning: the calo coverage behind the tracker is not complete! **" << endl; }
469
470 // final smeared EM energy // electromagnetic fraction F_em =1 for electrons;
471 float sE = currentCalo.getElectromagneticResolution().Smear(particle->E);
472 if (sE>0) {
473 sET = sE/cosh(particle->Eta);
474 // NB: ET is found via the calorimetry and not via the track curvature
475
476 TLorentzVector p;
477 p.SetPtEtaPhiE(sET,particle->Eta,particle->Phi,sE);
478 if (sET > DET->PTCUT_elec)
479 electron.push_back(D_Particle(p,particle->PID,particle->EtaCalo,particle->PhiCalo));
480 } else { sET=0;} // if negative smeared energy -- needed for the tracks
481 }
482 // 2.1a.2.3.4 Other charged particles : smear them for the tracks!
483 else { //other particles
484 D_CaloElement currentCalo = list_of_calorimeters.getElement(particle->EtaCalo);
485 float sEem = currentCalo.getElectromagneticResolution().Smear(particle->E * particle->getFem());
486 float sEhad = currentCalo.getHadronicResolution().Smear(particle->E * particle->getFhad());
487 float sE = ( (sEem>0)? sEem : 0 ) + ( (sEhad>0)? sEhad : 0 );
488 sET = sE/cosh(particle->EtaCalo);
489 }
490
491 // 2.1a.2.3.5 Tracks
492 if( (rand()%100) < DET->TRACK_eff && sET!=0) {
493 elementTrack = (TRootTracks*) branchTrack->NewEntry();
494 elementTrack->Set(particle->Eta, particle->Phi, particle->EtaCalo, particle->PhiCalo, sET, particle->Charge);
495 //TrackCentral.push_back(elementTrack->GetFourVector()); // tracks at vertex!
496 TrackCentral.push_back(*elementTrack); // tracks at vertex!
497 // TODO!!! apply a smearing on the position of the origin of the track
498 // TODO!!! elementTracks->SetPositionOut(Xout,Yout,Zout);
499 }
500 } // 2.1a.2.3 : if tracker/energy-flow
501 // 2.1a.2.4 Photons
502 // stays in the tracker -> track available -> gamma ID
503 else if( (pid == pGAMMA) && fabs(particle->EtaCalo)< DET->CEN_max_tracker ) {
504 float iEta=UNDEFINED, iPhi=UNDEFINED;
505 DET->BinEtaPhi(particle->PhiCalo,particle->EtaCalo,iPhi,iEta); // fills in iPhi and iEta
506 D_CaloTower tower(iEta,iPhi);
507 // stores the list of towers where to apply the photon ID algorithm. Just a trick for a faster search
508 list_of_towers_with_photon.addTower(tower);
509 }
510 // 2.1a.2.5 : very forward detectors
511 else if (DET->FLAG_vfd==1) {
512 // for the moment, only protons are transported
513 // BUT !!! could be a beam of other particles! (heavy ions?)
514 // BUT ALSO !!! if very forward muons, or others!
515 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
516 VFD->ZDC(treeWriter,branchZDC,particle);
517 }
518 // 2.1a.2.6: Zero degree calorimeter
519 //else if(DET->FLAG_zdc==1) {
520 //VFD->ZDC(treeWriter,branchZDC,particle);
521 //}
522
523 } // 2.1a.2 : if visible particle
524 delete particle;
525 } // loop on all particles 2.1a
526
527
528 // 2.1b loop on all (activated) towers
529 // at this stage, list_of_active_towers may contain several times the same tower
530 // first step is to merge identical towers, by matching their (iEta,iPhi)
531 list_of_active_towers.mergeDuplicates();
532 // Calotower smearing
533 list_of_active_towers.smearTowers(list_of_calorimeters);
534
535 for(unsigned int i=0; i<list_of_active_towers.size(); i++) {
536 float iEta = list_of_active_towers[i].getEta();
537 float iPhi = list_of_active_towers[i].getPhi();
538 float e = list_of_active_towers[i].getE();
539 if(iEta != UNDEFINED && iPhi != UNDEFINED && e!=0) {
540 elementCalo = (TRootCalo*) branchCalo->NewEntry();
541 elementCalo->set(list_of_active_towers[i]);
542 // not beautiful : should be improved!
543 TLorentzVector p;
544 p.SetPtEtaPhiE(list_of_active_towers[i].getET(), iEta, iPhi, e );
545 PhysicsTower Tower(LorentzVector(p.Px(),p.Py(),p.Pz(),p.E()));
546 towers.push_back(Tower);
547 }
548 } // loop on towers
549
550 // 2.1c photon ID
551 // list_of_towers_with_photon is the list of towers with photon candidates
552 // already smeared !
553 // sorts the vector and smears duplicates
554 list_of_towers_with_photon.mergeDuplicates();
555
556 for(unsigned int i=0; i<list_of_towers_with_photon.size(); i++) {
557 float eta = list_of_towers_with_photon[i].getEta();
558 float phi = list_of_towers_with_photon[i].getPhi();
559 D_CaloTower cal(list_of_active_towers.getElement(eta,phi));
560 if(cal.getEta() != UNDEFINED && cal.getPhi() != UNDEFINED && cal.getE() > 0) {
561 TLorentzVector p;
562 p.SetPtEtaPhiE(cal.getET(), eta,phi,cal.getE() );
563 if (cal.getET() > DET->PTCUT_gamma) { gamma.push_back(D_Particle(p,pGAMMA,p.Eta(),p.Phi())); }
564 }
565 } // for -- list of photons
566
567} // IF NEW ALGORITHM with energy flow
568
569
570
571
572 // 2.2 ********** Output preparation & complex objects ***********
573 // 2.2.1 ********************* sorting collections by decreasing pt
574 DET->SortedVector(electron);
575 for(unsigned int i=0; i < electron.size(); i++) {
576 elementElec = (TRootElectron*) branchElectron->NewEntry();
577 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
578 elementElec->EtaCalo = electron[i].EtaCalo();
579 elementElec->PhiCalo = electron[i].PhiCalo();
580 elementElec->Charge = sign(electron[i].PID());
581 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);//isolation based on tracks
582 } /////////////// HARDCODING
583 DET->SortedVector(muon);
584 for(unsigned int i=0; i < muon.size(); i++) {
585 elementMu = (TRootMuon*) branchMuon->NewEntry();
586 elementMu->Charge = sign(muon[i].PID());
587 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
588 elementMu->EtaCalo = muon[i].EtaCalo();
589 elementMu->PhiCalo = muon[i].PhiCalo();
590 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0); /////////////// HARDCODING
591 }
592 DET->SortedVector(gamma);
593 for(unsigned int i=0; i < gamma.size(); i++) {
594 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
595 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
596 }
597
598 // 2.2.2 ************* computes the Missing Transverse Momentum
599 TLorentzVector Att(0.,0.,0.,0.);
600 for(unsigned int i=0; i < towers.size(); i++)
601 {
602 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
603 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
604 {
605 PTmis = PTmis + Att;
606 // create a fastjet::PseudoJet with these components and put it onto
607 // back of the input_particles vector
608 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
609 }
610 }
611 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
612 elementEtmis->ET = (PTmis).Pt();
613 elementEtmis->Phi = (-PTmis).Phi();
614 elementEtmis->Px = (-PTmis).Px();
615 elementEtmis->Py = (-PTmis).Py();
616
617 // 2.2.3 ************* jets, B-tag, tau jets
618 sorted_jets=JETRUN->RunJets(input_particles);
619 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
620 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
621
622 treeWriter->Fill();
623 } // 2. Loop over all events ('for' loop)
624
625
626 cout <<"** **"<< endl;
627 cout <<"** Exiting detector simulation... **"<< endl;
628
629
630 treeWriter->Write();
631 delete treeWriter;
632 loopwatch.Stop();
633
634
635
636 // 3. ********** Trigger & Frog ***********
637 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
638 triggerwatch.Start();
639 if(DET->FLAG_trigger == 1)
640 {
641 cout <<"** **"<<endl;
642 cout <<"** ########### Start Trigger selection ########### **"<< endl;
643
644 // input
645 TChain chainT("Analysis");
646 chainT.Add(outputfilename.c_str());
647 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
648
649 // output
650 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
651 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
652 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
653 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
654 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
655 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
656
657 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
658 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
659
660
661 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
662 // loop on all entries
663 for(entryT = 0; entryT < allEntriesT; ++entryT) {
664 treeWriterT->Clear();
665 treeReaderT->ReadEntry(entryT);
666 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
667 treeWriterT->Fill();
668 } // loop on all entries
669 cout <<"** Exiting trigger simulation... **"<< endl;
670
671 treeWriterT->Write();
672 delete treeWriterT;
673 delete treeReaderT;
674 } // trigger
675 triggerwatch.Stop();
676
677
678 // 3.2 ************** FROG display
679 frogwatch.Start();
680 if(DET->FLAG_frog == 1) {
681 cout <<"** **"<<endl;
682 cout <<"** ################## Start FROG ################# **"<< endl;
683
684 FrogDisplay *FROG = new FrogDisplay(DET);
685 FROG->BuildEvents(outputfilename);
686 FROG->BuildGeom();
687 delete FROG;
688 }
689 frogwatch.Stop();
690
691
692
693
694// 4. ********** End & Exit ***********
695
696 globalwatch.Stop();
697 cout <<"** **"<< endl;
698 cout <<"** ################## Time report ################# **"<< endl;
699 cout << left << setw(32) <<"** Time report for "<<""
700 << left << setw(15) << allEntries <<""
701 << right << setw(22) <<"events **"<<endl;
702 cout <<"** **"<< endl;
703 cout << left << setw(10) <<"**"<<""
704 << left << setw(15) <<"Time (s):"<<""
705 << right << setw(15) <<"CPU"<<""
706 << right << setw(15) <<"Real"<<""
707 << right << setw(14) <<"**"<<endl;
708 cout << left << setw(10) <<"**"<<""
709 << left << setw(15) <<" + Global:"<<""
710 << right << setw(15) <<globalwatch.CpuTime()<<""
711 << right << setw(15) <<globalwatch.RealTime()<<""
712 << right << setw(14) <<"**"<<endl;
713 cout << left << setw(10) <<"**"<<""
714 << left << setw(15) <<" + Events:"<<""
715 << right << setw(15) <<loopwatch.CpuTime()<<""
716 << right << setw(15) <<loopwatch.RealTime()<<""
717 << right << setw(14) <<"**"<<endl;
718 if(DET->FLAG_trigger == 1)
719 {
720 cout << left << setw(10) <<"**"<<""
721 << left << setw(15) <<" + Trigger:"<<""
722 << right << setw(15) <<triggerwatch.CpuTime()<<""
723 << right << setw(15) <<triggerwatch.RealTime()<<""
724 << right << setw(14) <<"**"<<endl;
725 }
726 if(DET->FLAG_frog == 1)
727 {
728 cout << left << setw(10) <<"**"<<""
729 << left << setw(15) <<" + Frog:"<<""
730 << right << setw(15) <<frogwatch.CpuTime()<<""
731 << right << setw(15) <<frogwatch.RealTime()<<""
732 << right << setw(14) <<"**"<<endl;
733
734 }
735
736 cout <<"** **"<< endl;
737 cout <<"** Exiting Delphes ... **"<< endl;
738 cout <<"** **"<< endl;
739 cout <<"*********************************************************************"<< endl;
740 cout <<"*********************************************************************"<< endl;
741
742
743 delete treeReader;
744 delete DET;
745 delete TRIGT;
746 delete TRACP;
747 delete JETRUN;
748 delete VFD;
749 delete converter;
750
751// todo("TODO");
752}
Note: See TracBrowser for help on using the repository browser.