Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 243

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

include statements have been cleaned ; no explicit memory-leak

File size: 26.8 KB
RevLine 
[2]1/*
2 ---- Delphes ----
3 A Fast Simulator for general purpose LHC detector
4 S. Ovyn ~~~~ severine.ovyn@uclouvain.be
5
6 Center for Particle Physics and Phenomenology (CP3)
7 Universite Catholique de Louvain (UCL)
8 Louvain-la-Neuve, Belgium
9*/
10
11/// \file Delphes.cpp
12/// \brief executable for the Delphes
13
14#include "TChain.h"
15#include "TApplication.h"
[191]16#include "TStopwatch.h"
[228]17#include "TFile.h"
[2]18
[228]19#include "ExRootTreeReader.h"
20#include "ExRootTreeWriter.h"
21#include "ExRootTreeBranch.h"
[2]22
[228]23#include "DataConverter.h"
24#include "HEPEVTConverter.h"
25#include "LHEFConverter.h"
26#include "STDHEPConverter.h"
[2]27
[228]28#include "SmearUtil.h"
29#include "BFieldProp.h"
30#include "TriggerUtil.h"
31#include "VeryForward.h"
32#include "JetsUtil.h"
33#include "FrogUtil.h"
[2]34
[55]35#include <vector>
36#include <iostream>
[11]37
[2]38using namespace std;
39
40//------------------------------------------------------------------------------
41void todo(string filename) {
42 ifstream infile(filename.c_str());
43 cout << "** TODO list ..." << endl;
44 while(infile.good()) {
[94]45 string temp;
46 getline(infile,temp);
47 cout << "*" << temp << endl;
[2]48 }
49 cout << "** done...\n";
50}
51
52//------------------------------------------------------------------------------
53
54int main(int argc, char *argv[])
55{
56 int appargc = 2;
[228]57 char *appName= new char[20];
58 char *appOpt= new char[20];
59 sprintf(appName,"Delphes");
60 sprintf(appOpt,"-b");
61 char *appargv[] = {appName,appOpt};
[2]62 TApplication app(appName, &appargc, appargv);
[228]63 delete [] appName;
64 delete [] appOpt;
65
66
[71]67 if(argc != 4 && argc != 3 && argc != 5) {
[94]68 cout << " Usage: " << argv[0] << " input_file output_file [detector_card] [trigger_card] " << endl;
69 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
70 cout << " output_file - output file." << endl;
71 cout << " detector_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
72 cout << " trigger_card - Datacard containing the trigger algorithms (optional) "<<endl;
73 exit(1);
[2]74 }
[228]75
76
77
[212]78
79cout << endl << endl;
80
81cout <<"*********************************************************************"<< endl;
82cout <<"*********************************************************************"<< endl;
83cout <<"** **"<< endl;
84cout <<"** Welcome to **"<< endl;
85cout <<"** **"<< endl;
86cout <<"** **"<< endl;
87cout <<"** .ddddddd- lL hH **"<< endl;
88cout <<"** -Dd` `dD: Ll hH` **"<< endl;
89cout <<"** dDd dDd eeee. lL .pp+pp Hh+hhh` -eeee- `sssss **"<< endl;
90cout <<"** -Dd `DD ee. ee Ll .Pp. PP Hh. HH. ee. ee sSs **"<< endl;
91cout <<"** dD` dDd eEeee: lL. pP. pP hH hH` eEeee:` -sSSSs. **"<< endl;
92cout <<"** .Dd :dd eE. LlL PpppPP Hh Hh eE sSS **"<< endl;
93cout <<"** dddddd:. eee+: lL. pp. hh. hh eee+ sssssS **"<< endl;
94cout <<"** Pp **"<< endl;
95cout <<"** **"<< endl;
96cout <<"** Delphes, a framework for the fast simulation **"<< endl;
97cout <<"** of a generic collider experiment **"<< endl;
98cout <<"** **"<< endl;
99cout <<"** --- Version 1.3beta of Delphes --- **"<< endl;
100cout <<"** Last date of change: 29 January 2009 **"<< endl;
101cout <<"** **"<< endl;
102cout <<"** **"<< endl;
103cout <<"** This package uses: **"<< endl;
104cout <<"** ------------------ **"<< endl;
105cout <<"** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **"<< endl;
106cout <<"** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **"<< endl;
107cout <<"** FROG: **"<< endl;
108cout <<"** **"<< endl;
109cout <<"**-----------------------------------------------------------------**"<< endl;
110cout <<"** **"<< endl;
111cout <<"** Main authors: **"<< endl;
112cout <<"** ------------- **"<< endl;
113cout <<"** **"<< endl;
114cout <<"** Séverine Ovyn Xavier Rouby **"<< endl;
115cout <<"** severine.ovyn@uclouvain.be xavier.rouby@cern **"<< endl;
116cout <<"** Center for Particle Physics and Phenomenology (CP3) **"<< endl;
117cout <<"** Universite Catholique de Louvain (UCL) **"<< endl;
118cout <<"** Louvain-la-Neuve, Belgium **"<< endl;
119cout <<"** **"<< endl;
120cout <<"**-----------------------------------------------------------------**"<< endl;
121cout <<"** **"<< endl;
122cout <<"** Former Delphes versions and documentation can be found on : **"<< endl;
123cout <<"** http://www.fynu.ucl.ac.be/delphes.html **"<< endl;
124cout <<"** **"<< endl;
125cout <<"** **"<< endl;
126cout <<"** Disclaimer: this program is a beta version of Delphes and **"<< endl;
127cout <<"** therefore comes without guarantees. Beware of errors and please **"<< endl;
128cout <<"** give us your feedbacks about potential bugs **"<< endl;
129cout <<"** **"<< endl;
130cout <<"*********************************************************************"<< endl;
131cout <<"*********************************************************************"<< endl;
132
[178]133// 1. ********** initialisation ***********
134
[2]135 srand (time (NULL)); /* Initialisation du générateur */
[191]136 TStopwatch globalwatch, loopwatch, triggerwatch, frogwatch;
137 globalwatch.Start();
138
[2]139
140 //read the input TROOT file
141 string inputFileList(argv[1]), outputfilename(argv[2]);
142 if(outputfilename.find(".root") > outputfilename.length() ) {
[94]143 cout << "output_file should be a .root file!\n";
144 exit(1);
[2]145 }
[44]146 //create output log-file name
[45]147 string forLog = outputfilename;
148 string LogName = forLog.erase(forLog.find(".root"));
[44]149 LogName = LogName+"_run.log";
[94]150
[2]151 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
152 outputFile->Close();
[94]153
[2]154 string line;
155 ifstream infile(inputFileList.c_str());
156 infile >> line; // the first line determines the type of input files
[94]157
[44]158 //read the datacard input file
159 string DetDatacard("");
[212]160 if(argc>=4)DetDatacard =argv[3];
[94]161
[55]162 //Smearing information
[44]163 RESOLution *DET = new RESOLution();
[212]164 cout <<"** **"<< endl;
165 cout <<"** ####### Start reading DETECTOR parameters ####### **"<< endl;
166 cout << left << setw(40) <<"** Opening configuration card: "<<""
167 << left << setw(20) << DetDatacard <<""
168 << right << setw(9) <<" **"<<""<<endl;
[44]169 DET->ReadDataCard(DetDatacard);
170 DET->Logfile(LogName);
[212]171 cout << left << setw(42) <<"** Parameters summarised in file: "<<""
172 << left << setw(20) << LogName <<""
173 << right << setw(7) <<"**"<<""<<endl;
174 cout <<"** **"<< endl;
[228]175
[79]176 //read the trigger input file
[80]177 string TrigDatacard("data/trigger.dat");
[79]178 if(argc==5) TrigDatacard =argv[4];
[94]179
[55]180 //Trigger information
[72]181 TriggerTable *TRIGT = new TriggerTable();
[80]182 TRIGT->TriggerCardReader(TrigDatacard.c_str());
[72]183 TRIGT->PrintTriggerTable(LogName);
[212]184 if(DET->FLAG_trigger == 1)
185 {
186 cout <<"** ########### Start reading TRIGGER card ########## **"<< endl;
187 cout << left << setw(40) <<"** Opening configuration card: "<<""
188 << left << setw(20) << TrigDatacard <<""
189 << right << setw(9) <<" **"<<""<<endl;
190 cout <<"** **"<< endl;
191 }
[94]192
[55]193 //Propagation of tracks in the B field
[228]194 TrackPropagation *TRACP = new TrackPropagation(DET);
195 //TrackPropagation *TRACP = new TrackPropagation(DetDatacard);
[94]196
[55]197 //Jet information
[228]198 JetsUtil *JETRUN = new JetsUtil(DET);
199 //JetsUtil *JETRUN = new JetsUtil(DetDatacard);
[94]200
[55]201 //VFD information
[228]202 VeryForward * VFD = new VeryForward(DET);
203 //VeryForward * VFD = new VeryForward(DetDatacard);
[178]204
205 // data converters
[228]206 DataConverter *converter=NULL;
[212]207 cout <<"** **"<<endl;
208 cout <<"** ####### Start convertion to TRoot format ######## **"<< endl;
209
[2]210 if(strstr(line.c_str(),".hep"))
211 {
[212]212 cout <<"** StdHEP file format detected **"<<endl;
213 cout <<"** This can take several minutes **"<< endl;
[2]214 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
215 }
216 else if(strstr(line.c_str(),".lhe"))
217 {
[212]218 cout <<"** LHEF file format detected **"<<endl;
219 cout <<"** This can take several minutes **"<< endl;
[2]220 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
221 }
222 else if(strstr(line.c_str(),".root"))
223 {
[212]224 cout <<"** h2root file format detected **"<<endl;
225 cout <<"** This can take several minutes **"<< endl;
[2]226 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
227 }
[212]228 else {
229 cout << left << setw(4) <<"** "<<""
230 << left << setw(63) << line.c_str() <<""
231 << right << setw(2) <<"**"<<endl;
232 cout <<"** File format not identified -- Exiting... **"<< endl;
233 cout <<"** **"<< endl;
234 cout <<"*********************************************************************"<< endl;
235 return -1;};
236 cout <<"** Exiting conversion... **"<< endl;
237
[2]238 TChain chain("GEN");
239 chain.Add(outputfilename.c_str());
240 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
241 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
242 TIter itGen((TCollection*)branchGen);
243
[178]244 //Output file : contents of the analysis object data
[2]245 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
246 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
247 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
248 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
249 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
250 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
251 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
252 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
253 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
254 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
255 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
256 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
[30]257
[2]258 TRootGenParticle *particle;
259 TRootETmis *elementEtmis;
260 TRootElectron *elementElec;
261 TRootMuon *elementMu;
262 TRootPhoton *elementPhoton;
263 TRootTracks *elementTracks;
264 TRootCalo *elementCalo;
265
[184]266 TLorentzVector genMomentum(0,0,0,0); // four-momentum at the vertex
267 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
268 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
[2]269 LorentzVector jetMomentum;
[94]270
[55]271 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
272 vector<fastjet::PseudoJet> sorted_jets;
[2]273 vector<TLorentzVector> TrackCentral;
274 vector<PhysicsTower> towers;
[73]275 vector<ParticleUtil> electron;
276 vector<ParticleUtil> muon;
[74]277 vector<ParticleUtil> gamma;
[98]278
[30]279 TSimpleArray<TRootGenParticle> NFCentralQ;
[100]280 float iPhi=0,iEta=0;
281
[178]282
283
284// 2. ********** Loop over all events ***********
[2]285 Long64_t entry, allEntries = treeReader->GetEntries();
[212]286 cout <<"** **"<<endl;
287 cout <<"** ####### Start fast detector simulation ######## **"<< endl;
288 cout << left << setw(52) <<"** Total number of events to run: "<<""
289 << left << setw(15) << allEntries <<""
290 << right << setw(2) <<"**"<<endl;
291
[191]292 loopwatch.Start();
[178]293
294 // loop on all events
295 for(entry = 0; entry < allEntries; ++entry)
[2]296 {
297 TLorentzVector PTmis(0,0,0,0);
298 treeReader->ReadEntry(entry);
299 treeWriter->Clear();
[228]300 if((entry % 100) == 0 && entry > 0 ) {
[212]301 cout << left << setw(52) <<"** Processing element # "<<""
302 << left << setw(15) << entry <<""
303 << right << setw(2) <<"**"<<endl;
304 }
[2]305
[30]306 electron.clear();
307 muon.clear();
[74]308 gamma.clear();
[30]309 NFCentralQ.Clear();
310
[2]311 TrackCentral.clear();
312 towers.clear();
[11]313 input_particles.clear();
[30]314
[178]315 // 2.1 Loop over all particles in event
[74]316 itGen.Reset();
[178]317 while( (particle = (TRootGenParticle*) itGen.Next()) ) {
[55]318 int pid = abs(particle->PID);
[178]319
320
321 // 2.1.1********************* preparation for the b-tagging
[94]322 //// This subarray is needed for the B-jet algorithm
[2]323 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
324 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 ?
[94]325 fabs(particle->Eta) < DET->CEN_max_tracker &&
[2]326 particle->Status != 1 &&
327 particle->PT > DET->PT_QUARKS_MIN ) {
328 NFCentralQ.Add(particle);
329 }
[94]330
[178]331 // 2.1.2 ********************* central detector: keep only visible particles
[2]332 // keeps only final particles, visible by the central detector, including the fiducial volume
333 // the ordering of conditions have been optimised for speed : put first the STATUS condition
334 if( (particle->Status == 1) &&
[59]335 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
[94]336 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
337 )
[59]338 {
[94]339 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
[184]340 genMomentumBfield = genMomentum;
[178]341
342 // 2.1.2.1 ********************* central detector: magnetic field
[184]343 // genMomentumBfield is then changed with respect to the magnetic field
344 if(DET->FLAG_bfield==1) TRACP->Propagation(particle,genMomentumBfield);
345 float eta=fabs(genMomentumBfield.Eta());
[178]346
347
348 // 2.1.2.2 ********************* central detector: smearing (calorimeters, muon chambers)
[94]349 switch(pid) {
[74]350
[94]351 case pE: // all electrons with eta < DET->MAX_CALO_FWD
[184]352 DET->SmearElectron(genMomentumBfield);
353 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_elec){
354 electron.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]355 }
356 break; // case pE
357 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
[184]358 DET->SmearElectron(genMomentumBfield);
359 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_gamma) {
360 gamma.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]361 }
362 break; // case pGAMMA
363 case pMU: // all muons with eta < DET->MAX_MU
[184]364 DET->SmearMu(genMomentumBfield);
365 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_mu && genMomentumBfield.Pt() > DET->PTCUT_muon){
366 muon.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]367 }
368 break; // case pMU
369 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
370 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
[184]371 DET->SmearHadron(genMomentumBfield, 0.7);
[94]372 break; // case hadron
373 default: // all other final particles with eta < DET->MAX_CALO_FWD
[184]374 DET->SmearHadron(genMomentumBfield, 1.0);
[94]375 break;
[178]376 } // 2.1.2.2 switch (pid)
377
378
379 // 2.1.2.3 ********************* central detector: calotowers
[94]380 // all final particles but muons and neutrinos
[100]381 // for calorimetric towers and missing PT
[94]382 int charge=Charge(pid);
[184]383 if(genMomentumBfield.E() !=0 && pid != pMU) {
[178]384 // in case the Bfield is not simulated, checks that charged particles have enough pt to reach the calos
[184]385 if ( !DET->FLAG_bfield && charge!=0 && genMomentumBfield.Pt() <= DET->TRACK_ptmin ) { /* particules do not reach calos */ }
[178]386 else { // particles reach calos
387 // applies the calo segmentation and returns iEta & iPhi
[184]388 DET->BinEtaPhi(genMomentumBfield.Phi(), genMomentumBfield.Eta(), iPhi, iEta);
[178]389 if(iEta != -100 && iPhi != -100) {
[184]390 momentumCaloSegmentation.SetPtEtaPhiE(genMomentumBfield.Pt(),iEta,iPhi,genMomentumBfield.E());
[100]391 elementCalo = (TRootCalo*) branchCalo->NewEntry();
[184]392 elementCalo->Set(momentumCaloSegmentation);
393 PhysicsTower Tower(LorentzVector(momentumCaloSegmentation.Px(),momentumCaloSegmentation.Py(),momentumCaloSegmentation.Pz(),momentumCaloSegmentation.E()));
[100]394 towers.push_back(Tower);
[178]395 } // if iEta != -100
396 } // else : when particles reach the calos
397 } // 2.1.2.3 calotowers
[94]398
[178]399
400 // 2.1.2.4 ********************* central detector: tracks
[94]401 // all final charged particles
402 if(
[184]403 (genMomentumBfield.E()!=0) &&
404 (fabs(genMomentumBfield.Eta()) < DET->CEN_max_tracker) &&
405 (DET->FLAG_bfield || ( !DET->FLAG_bfield && genMomentumBfield.Pt() > DET->TRACK_ptmin )) &&
[178]406 // if bfield not simulated, pt should be high enough to be taken into account
[94]407 ((rand()%100) < DET->TRACK_eff) &&
408 (charge!=0)
409 )
410 {
411 elementTracks = (TRootTracks*) branchTracks->NewEntry();
[184]412 elementTracks->Set(genMomentum); // fills px,py,pz,pt,e,eta,phi at vertex
413 elementTracks->Etaout = genMomentumBfield.Eta();
414 elementTracks->Phiout = genMomentumBfield.Phi();
415 // TODO!!! apply a smearing on the position of the origin of the track
416 // elementTracks->SetPosition(particle->X,particle->Y,particle->Z);
417 // uses the output of the bfield computation : Xout=... Yout=... Zout...
418 // TODO!!! elementTrakcs->SetPositionOut(Xout,Yout,Zout);
419 TrackCentral.push_back(genMomentum); // tracks at vertex!
[178]420 } // 2.1.2.4 tracks
[94]421
[178]422 } // 2.1.2 central detector
[30]423
[178]424 // 2.1.3 ********************* forward detectors: zdc
425 if(DET->FLAG_vfd==1) {
[100]426 VFD->ZDC(treeWriter,branchZDC,particle);
427 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
[178]428 }
[11]429
[178]430 } // 2.1 while : loop on all particles of the event.
431
432
433
434 // 2.2 ********** Output preparation & complex objects ***********
435
436 // 2.2.1 ********************* sorting collections by decreasing pt
[74]437 DET->SortedVector(electron);
438 for(unsigned int i=0; i < electron.size(); i++) {
439 elementElec = (TRootElectron*) branchElectron->NewEntry();
440 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
441 elementElec->Charge = sign(electron[i].PID());
[184]442 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);//isolation based on tracks
[30]443 }
[74]444 DET->SortedVector(muon);
[30]445 for(unsigned int i=0; i < muon.size(); i++) {
[74]446 elementMu = (TRootMuon*) branchMuon->NewEntry();
447 elementMu->Charge = sign(muon[i].PID());
448 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
449 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
[30]450 }
[74]451 DET->SortedVector(gamma);
452 for(unsigned int i=0; i < gamma.size(); i++) {
453 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
454 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
455 }
[30]456
[178]457 // 2.2.2 ************* computes the Missing Transverse Momentum
[71]458 TLorentzVector Att(0.,0.,0.,0.);
459 for(unsigned int i=0; i < towers.size(); i++)
460 {
[107]461 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
462 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
463 {
464 PTmis = PTmis + Att;
465 // create a fastjet::PseudoJet with these components and put it onto
466 // back of the input_particles vector
467 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
468 }
[71]469 }
470 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
471 elementEtmis->ET = (PTmis).Pt();
472 elementEtmis->Phi = (-PTmis).Phi();
473 elementEtmis->Px = (-PTmis).Px();
474 elementEtmis->Py = (-PTmis).Py();
[74]475
[178]476 // 2.2.3 ************* B-tag, tau jets
[55]477 sorted_jets=JETRUN->RunJets(input_particles);
478 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
479 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
[74]480
[72]481 treeWriter->Fill();
[178]482 } // 2. Loop over all events ('for' loop)
[212]483
484 cout <<"** **"<< endl;
485 cout <<"** Exiting detector simulation... **"<< endl;
486
487
[2]488 treeWriter->Write();
[77]489 delete treeWriter;
[191]490 loopwatch.Stop();
[178]491
492
493
[212]494 // 3. ********** Trigger & Frog ***********
[178]495 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
[191]496 triggerwatch.Start();
[94]497 if(DET->FLAG_trigger == 1)
[72]498 {
[212]499 cout <<"** **"<<endl;
500 cout <<"** ########### Start Trigger selection ########### **"<< endl;
501
[178]502 // input
[72]503 TChain chainT("Analysis");
504 chainT.Add(outputfilename.c_str());
505 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
[74]506
[178]507 // output
[72]508 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
509 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
510 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
511 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
512 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
513 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
[74]514
[72]515 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
516 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
[178]517
518
[72]519 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
[178]520 // loop on all entries
521 for(entryT = 0; entryT < allEntriesT; ++entryT) {
[72]522 treeWriterT->Clear();
523 treeReaderT->ReadEntry(entryT);
524 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
525 treeWriterT->Fill();
[178]526 } // loop on all entries
[212]527 cout <<"** Exiting trigger simulation... **"<< endl;
[94]528
[72]529 treeWriterT->Write();
530 delete treeWriterT;
[228]531 delete treeReaderT;
[178]532 } // trigger
[191]533 triggerwatch.Stop();
[178]534
535
536 // 3.2 ************** FROG display
[191]537 frogwatch.Start();
[228]538 if(DET->FLAG_frog == 1) {
[212]539 cout <<"** **"<<endl;
540 cout <<"** ################## Start FROG ################# **"<< endl;
541
[228]542 FrogDisplay *FROG = new FrogDisplay(DET);
543 FROG->BuildEvents(outputfilename);
544 FROG->BuildGeom();
545 delete FROG;
546 }
547 frogwatch.Stop();
[94]548
[178]549
550
551
552// 4. ********** End & Exit ***********
[212]553
[191]554 globalwatch.Stop();
[212]555 cout <<"** **"<< endl;
556 cout <<"** ################## Time report ################# **"<< endl;
557 cout << left << setw(32) <<"** Time report for "<<""
558 << left << setw(15) << allEntries <<""
559 << right << setw(22) <<"events **"<<endl;
560 cout <<"** **"<< endl;
561 cout << left << setw(10) <<"**"<<""
562 << left << setw(15) <<"Time (s):"<<""
563 << right << setw(15) <<"CPU"<<""
564 << right << setw(15) <<"Real"<<""
565 << right << setw(14) <<"**"<<endl;
566 cout << left << setw(10) <<"**"<<""
567 << left << setw(15) <<" + Global:"<<""
568 << right << setw(15) <<globalwatch.CpuTime()<<""
569 << right << setw(15) <<globalwatch.RealTime()<<""
570 << right << setw(14) <<"**"<<endl;
571 cout << left << setw(10) <<"**"<<""
572 << left << setw(15) <<" + Events:"<<""
573 << right << setw(15) <<loopwatch.CpuTime()<<""
574 << right << setw(15) <<loopwatch.RealTime()<<""
575 << right << setw(14) <<"**"<<endl;
[191]576 if(DET->FLAG_trigger == 1)
[212]577 {
578 cout << left << setw(10) <<"**"<<""
579 << left << setw(15) <<" + Trigger:"<<""
580 << right << setw(15) <<triggerwatch.CpuTime()<<""
581 << right << setw(15) <<triggerwatch.RealTime()<<""
582 << right << setw(14) <<"**"<<endl;
583 }
[191]584 if(DET->FLAG_frog == 1)
[212]585 {
586 cout << left << setw(10) <<"**"<<""
587 << left << setw(15) <<" + Frog:"<<""
588 << right << setw(15) <<frogwatch.CpuTime()<<""
589 << right << setw(15) <<frogwatch.RealTime()<<""
590 << right << setw(14) <<"**"<<endl;
[191]591
[212]592 }
[2]593
[212]594 cout <<"** **"<< endl;
595 cout <<"** Exiting Delphes ... **"<< endl;
596 cout <<"** **"<< endl;
597 cout <<"*********************************************************************"<< endl;
598 cout <<"*********************************************************************"<< endl;
599
600
[2]601 delete treeReader;
602 delete DET;
[74]603 delete TRIGT;
604 delete TRACP;
605 delete JETRUN;
606 delete VFD;
[228]607 delete converter;
[94]608
[191]609// todo("TODO");
[2]610}
Note: See TracBrowser for help on using the repository browser.