Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 264

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

first test 2.0

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