Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 194

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

Performance of timing included

File size: 18.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"
[2]17
18#include "Utilities/ExRootAnalysis/interface/ExRootTreeReader.h"
19#include "Utilities/ExRootAnalysis/interface/ExRootTreeWriter.h"
20#include "Utilities/ExRootAnalysis/interface/ExRootTreeBranch.h"
21
22#include "interface/DataConverter.h"
23#include "interface/HEPEVTConverter.h"
24#include "interface/LHEFConverter.h"
25#include "interface/STDHEPConverter.h"
26
27#include "interface/SmearUtil.h"
[55]28#include "interface/BFieldProp.h"
29#include "interface/TriggerUtil.h"
30#include "interface/VeryForward.h"
31#include "interface/JetUtils.h"
[94]32#include "interface/FrogUtil.h"
[2]33
[55]34#include <vector>
35#include <iostream>
[11]36
[2]37using namespace std;
38
39//------------------------------------------------------------------------------
40void todo(string filename) {
41 ifstream infile(filename.c_str());
42 cout << "** TODO list ..." << endl;
43 while(infile.good()) {
[94]44 string temp;
45 getline(infile,temp);
46 cout << "*" << temp << endl;
[2]47 }
48 cout << "** done...\n";
49}
50
51//------------------------------------------------------------------------------
52
53int main(int argc, char *argv[])
54{
55 int appargc = 2;
[55]56 char *appName = "Delphes";
[2]57 char *appargv[] = {appName, "-b"};
58 TApplication app(appName, &appargc, appargv);
[94]59
[71]60 if(argc != 4 && argc != 3 && argc != 5) {
[94]61 cout << " Usage: " << argv[0] << " input_file output_file [detector_card] [trigger_card] " << endl;
62 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
63 cout << " output_file - output file." << endl;
64 cout << " detector_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
65 cout << " trigger_card - Datacard containing the trigger algorithms (optional) "<<endl;
66 exit(1);
[2]67 }
[94]68
[178]69// 1. ********** initialisation ***********
70
[2]71 srand (time (NULL)); /* Initialisation du générateur */
[191]72 TStopwatch globalwatch, loopwatch, triggerwatch, frogwatch;
73 globalwatch.Start();
74
[2]75
76 //read the input TROOT file
77 string inputFileList(argv[1]), outputfilename(argv[2]);
78 if(outputfilename.find(".root") > outputfilename.length() ) {
[94]79 cout << "output_file should be a .root file!\n";
80 exit(1);
[2]81 }
[44]82 //create output log-file name
[45]83 string forLog = outputfilename;
84 string LogName = forLog.erase(forLog.find(".root"));
[44]85 LogName = LogName+"_run.log";
[94]86
[2]87 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
88 outputFile->Close();
[94]89
[2]90 string line;
91 ifstream infile(inputFileList.c_str());
92 infile >> line; // the first line determines the type of input files
[94]93
[44]94 //read the datacard input file
95 string DetDatacard("");
[94]96 if(argc>=4) DetDatacard =argv[3];
97
[55]98 //Smearing information
[44]99 RESOLution *DET = new RESOLution();
100 DET->ReadDataCard(DetDatacard);
101 DET->Logfile(LogName);
[94]102
[79]103 //read the trigger input file
[80]104 string TrigDatacard("data/trigger.dat");
[79]105 if(argc==5) TrigDatacard =argv[4];
[94]106
[55]107 //Trigger information
[72]108 TriggerTable *TRIGT = new TriggerTable();
[80]109 TRIGT->TriggerCardReader(TrigDatacard.c_str());
[72]110 TRIGT->PrintTriggerTable(LogName);
[94]111
[55]112 //Propagation of tracks in the B field
[100]113 TrackPropagation *TRACP = new TrackPropagation(DetDatacard);
[94]114
[55]115 //Jet information
[100]116 JetsUtil *JETRUN = new JetsUtil(DetDatacard);
[94]117
[55]118 //VFD information
[100]119 VeryForward * VFD = new VeryForward(DetDatacard);
[178]120
121 // data converters
[2]122 DataConverter *converter=0;
[94]123
[2]124 if(strstr(line.c_str(),".hep"))
125 {
[44]126 cout<<"#**********************************************************************"<<endl;
127 cout<<"#********** StdHEP file format detected *************"<<endl;
128 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
129 cout<<"#**********************************************************************"<<endl;
[2]130 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
131 }
132 else if(strstr(line.c_str(),".lhe"))
133 {
[44]134 cout<<"#**********************************************************************"<<endl;
135 cout<<"#*********** LHEF file format detected ************"<<endl;
136 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
137 cout<<"#**********************************************************************"<<endl;
[2]138 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
139 }
140 else if(strstr(line.c_str(),".root"))
141 {
[44]142 cout<<"#**********************************************************************"<<endl;
143 cout<<"#********** h2root file format detected *************"<<endl;
144 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
145 cout<<"#**********************************************************************"<<endl;
[2]146 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
147 }
148 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
[30]149
[2]150 TChain chain("GEN");
151 chain.Add(outputfilename.c_str());
152 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
153 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
154 TIter itGen((TCollection*)branchGen);
155
[178]156 //Output file : contents of the analysis object data
[2]157 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
158 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
159 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
160 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
161 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
162 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
163 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
164 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
165 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
166 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
167 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
168 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
[30]169
[2]170 TRootGenParticle *particle;
171 TRootETmis *elementEtmis;
172 TRootElectron *elementElec;
173 TRootMuon *elementMu;
174 TRootPhoton *elementPhoton;
175 TRootTracks *elementTracks;
176 TRootCalo *elementCalo;
177
[184]178 TLorentzVector genMomentum(0,0,0,0); // four-momentum at the vertex
179 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
180 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
[2]181 LorentzVector jetMomentum;
[94]182
[55]183 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
184 vector<fastjet::PseudoJet> sorted_jets;
[2]185 vector<TLorentzVector> TrackCentral;
186 vector<PhysicsTower> towers;
[73]187 vector<ParticleUtil> electron;
188 vector<ParticleUtil> muon;
[74]189 vector<ParticleUtil> gamma;
[98]190
[30]191 TSimpleArray<TRootGenParticle> NFCentralQ;
[100]192 float iPhi=0,iEta=0;
193
[178]194
195
196// 2. ********** Loop over all events ***********
[2]197 Long64_t entry, allEntries = treeReader->GetEntries();
[178]198 cout << "** The input list contains " << allEntries << " events" << endl;
[191]199 loopwatch.Start();
[178]200
201 // loop on all events
202 for(entry = 0; entry < allEntries; ++entry)
[2]203 {
204 TLorentzVector PTmis(0,0,0,0);
205 treeReader->ReadEntry(entry);
206 treeWriter->Clear();
207 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
208
[30]209 electron.clear();
210 muon.clear();
[74]211 gamma.clear();
[30]212 NFCentralQ.Clear();
213
[2]214 TrackCentral.clear();
215 towers.clear();
[11]216 input_particles.clear();
[30]217
[178]218 // 2.1 Loop over all particles in event
[74]219 itGen.Reset();
[178]220 while( (particle = (TRootGenParticle*) itGen.Next()) ) {
[55]221 int pid = abs(particle->PID);
[178]222
223
224 // 2.1.1********************* preparation for the b-tagging
[94]225 //// This subarray is needed for the B-jet algorithm
[2]226 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
227 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]228 fabs(particle->Eta) < DET->CEN_max_tracker &&
[2]229 particle->Status != 1 &&
230 particle->PT > DET->PT_QUARKS_MIN ) {
231 NFCentralQ.Add(particle);
232 }
[94]233
[178]234 // 2.1.2 ********************* central detector: keep only visible particles
[2]235 // keeps only final particles, visible by the central detector, including the fiducial volume
236 // the ordering of conditions have been optimised for speed : put first the STATUS condition
237 if( (particle->Status == 1) &&
[59]238 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
[94]239 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
240 )
[59]241 {
[94]242 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
[184]243 genMomentumBfield = genMomentum;
[178]244
245 // 2.1.2.1 ********************* central detector: magnetic field
[184]246 // genMomentumBfield is then changed with respect to the magnetic field
247 if(DET->FLAG_bfield==1) TRACP->Propagation(particle,genMomentumBfield);
248 float eta=fabs(genMomentumBfield.Eta());
[178]249
250
251 // 2.1.2.2 ********************* central detector: smearing (calorimeters, muon chambers)
[94]252 switch(pid) {
[74]253
[94]254 case pE: // all electrons with eta < DET->MAX_CALO_FWD
[184]255 DET->SmearElectron(genMomentumBfield);
256 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_elec){
257 electron.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]258 }
259 break; // case pE
260 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
[184]261 DET->SmearElectron(genMomentumBfield);
262 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_gamma) {
263 gamma.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]264 }
265 break; // case pGAMMA
266 case pMU: // all muons with eta < DET->MAX_MU
[184]267 DET->SmearMu(genMomentumBfield);
268 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_mu && genMomentumBfield.Pt() > DET->PTCUT_muon){
269 muon.push_back(ParticleUtil(genMomentumBfield,particle->PID));
[94]270 }
271 break; // case pMU
272 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
273 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
[184]274 DET->SmearHadron(genMomentumBfield, 0.7);
[94]275 break; // case hadron
276 default: // all other final particles with eta < DET->MAX_CALO_FWD
[184]277 DET->SmearHadron(genMomentumBfield, 1.0);
[94]278 break;
[178]279 } // 2.1.2.2 switch (pid)
280
281
282 // 2.1.2.3 ********************* central detector: calotowers
[94]283 // all final particles but muons and neutrinos
[100]284 // for calorimetric towers and missing PT
[94]285 int charge=Charge(pid);
[184]286 if(genMomentumBfield.E() !=0 && pid != pMU) {
[178]287 // in case the Bfield is not simulated, checks that charged particles have enough pt to reach the calos
[184]288 if ( !DET->FLAG_bfield && charge!=0 && genMomentumBfield.Pt() <= DET->TRACK_ptmin ) { /* particules do not reach calos */ }
[178]289 else { // particles reach calos
290 // applies the calo segmentation and returns iEta & iPhi
[184]291 DET->BinEtaPhi(genMomentumBfield.Phi(), genMomentumBfield.Eta(), iPhi, iEta);
[178]292 if(iEta != -100 && iPhi != -100) {
[184]293 momentumCaloSegmentation.SetPtEtaPhiE(genMomentumBfield.Pt(),iEta,iPhi,genMomentumBfield.E());
[100]294 elementCalo = (TRootCalo*) branchCalo->NewEntry();
[184]295 elementCalo->Set(momentumCaloSegmentation);
296 PhysicsTower Tower(LorentzVector(momentumCaloSegmentation.Px(),momentumCaloSegmentation.Py(),momentumCaloSegmentation.Pz(),momentumCaloSegmentation.E()));
[100]297 towers.push_back(Tower);
[178]298 } // if iEta != -100
299 } // else : when particles reach the calos
300 } // 2.1.2.3 calotowers
[94]301
[178]302
303 // 2.1.2.4 ********************* central detector: tracks
[94]304 // all final charged particles
305 if(
[184]306 (genMomentumBfield.E()!=0) &&
307 (fabs(genMomentumBfield.Eta()) < DET->CEN_max_tracker) &&
308 (DET->FLAG_bfield || ( !DET->FLAG_bfield && genMomentumBfield.Pt() > DET->TRACK_ptmin )) &&
[178]309 // if bfield not simulated, pt should be high enough to be taken into account
[94]310 ((rand()%100) < DET->TRACK_eff) &&
311 (charge!=0)
312 )
313 {
314 elementTracks = (TRootTracks*) branchTracks->NewEntry();
[184]315 elementTracks->Set(genMomentum); // fills px,py,pz,pt,e,eta,phi at vertex
316 elementTracks->Etaout = genMomentumBfield.Eta();
317 elementTracks->Phiout = genMomentumBfield.Phi();
318 // TODO!!! apply a smearing on the position of the origin of the track
319 // elementTracks->SetPosition(particle->X,particle->Y,particle->Z);
320 // uses the output of the bfield computation : Xout=... Yout=... Zout...
321 // TODO!!! elementTrakcs->SetPositionOut(Xout,Yout,Zout);
322 TrackCentral.push_back(genMomentum); // tracks at vertex!
[178]323 } // 2.1.2.4 tracks
[94]324
[178]325 } // 2.1.2 central detector
[30]326
[178]327 // 2.1.3 ********************* forward detectors: zdc
328 if(DET->FLAG_vfd==1) {
[100]329 VFD->ZDC(treeWriter,branchZDC,particle);
330 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
[178]331 }
[11]332
[178]333 } // 2.1 while : loop on all particles of the event.
334
335
336
337 // 2.2 ********** Output preparation & complex objects ***********
338
339 // 2.2.1 ********************* sorting collections by decreasing pt
[74]340 DET->SortedVector(electron);
341 for(unsigned int i=0; i < electron.size(); i++) {
342 elementElec = (TRootElectron*) branchElectron->NewEntry();
343 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
344 elementElec->Charge = sign(electron[i].PID());
[184]345 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);//isolation based on tracks
[30]346 }
[74]347 DET->SortedVector(muon);
[30]348 for(unsigned int i=0; i < muon.size(); i++) {
[74]349 elementMu = (TRootMuon*) branchMuon->NewEntry();
350 elementMu->Charge = sign(muon[i].PID());
351 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
352 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
[30]353 }
[74]354 DET->SortedVector(gamma);
355 for(unsigned int i=0; i < gamma.size(); i++) {
356 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
357 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
358 }
[30]359
[178]360 // 2.2.2 ************* computes the Missing Transverse Momentum
[71]361 TLorentzVector Att(0.,0.,0.,0.);
362 for(unsigned int i=0; i < towers.size(); i++)
363 {
[107]364 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
365 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
366 {
367 PTmis = PTmis + Att;
368 // create a fastjet::PseudoJet with these components and put it onto
369 // back of the input_particles vector
370 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
371 }
[71]372 }
373 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
374 elementEtmis->ET = (PTmis).Pt();
375 elementEtmis->Phi = (-PTmis).Phi();
376 elementEtmis->Px = (-PTmis).Px();
377 elementEtmis->Py = (-PTmis).Py();
[74]378
[178]379 // 2.2.3 ************* B-tag, tau jets
[55]380 sorted_jets=JETRUN->RunJets(input_particles);
381 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
382 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
[74]383
[72]384 treeWriter->Fill();
[178]385 } // 2. Loop over all events ('for' loop)
[74]386
[2]387 treeWriter->Write();
[77]388 delete treeWriter;
[191]389 loopwatch.Stop();
[178]390
391
392
393// 3. ********** Trigger & Frog ***********
394 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
[191]395 triggerwatch.Start();
[94]396 if(DET->FLAG_trigger == 1)
[72]397 {
[178]398 // input
[72]399 TChain chainT("Analysis");
400 chainT.Add(outputfilename.c_str());
401 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
[74]402
[178]403 // output
[72]404 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
405 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
406 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
407 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
408 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
409 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
[74]410
[72]411 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
412 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
[178]413
414
[72]415 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
[178]416 cout << "** Trigger: the 'Analysis' tree contains " << allEntriesT << " events" << endl;
417 // loop on all entries
418 for(entryT = 0; entryT < allEntriesT; ++entryT) {
[72]419 treeWriterT->Clear();
420 treeReaderT->ReadEntry(entryT);
421 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
422 treeWriterT->Fill();
[178]423 } // loop on all entries
[94]424
[72]425 treeWriterT->Write();
426 delete treeWriterT;
[178]427 } // trigger
[191]428 triggerwatch.Stop();
[178]429
430
431 // 3.2 ************** FROG display
[191]432 frogwatch.Start();
[94]433 if(DET->FLAG_frog == 1)
434 {
435 FrogDisplay *FROG = new FrogDisplay();
436 FROG->BuidEvents(outputfilename,DET->NEvents_Frog);
[98]437 FROG->BuildGeom(DetDatacard);
[94]438 }
[191]439 frogwatch.Stop();
[94]440
[178]441
442
443
444// 4. ********** End & Exit ***********
[2]445 cout << "** Exiting..." << endl;
[191]446 globalwatch.Stop();
447 cout << "** Time report for " << allEntries << " events.\n";
448 cout << " + Time (s): \tCPU \t real"<< endl;
449 cout << " + Global: \t" << globalwatch.CpuTime() << " \t " << globalwatch.RealTime() << endl;
450 cout << " + Events: \t" << loopwatch.CpuTime() << " \t " << loopwatch.RealTime() << endl;
451 if(DET->FLAG_trigger == 1)
452 cout << " + Trigger: \t" << triggerwatch.CpuTime() << " \t " << triggerwatch.RealTime() << endl;
453 if(DET->FLAG_frog == 1)
454 cout << " + Frog: \t" << frogwatch.CpuTime() << " \t " << frogwatch.RealTime() << endl;
455
[2]456
457 delete treeReader;
458 delete DET;
[74]459 delete TRIGT;
460 delete TRACP;
461 delete JETRUN;
462 delete VFD;
[2]463 if(converter) delete converter;
[94]464
[191]465// todo("TODO");
[2]466}
Note: See TracBrowser for help on using the repository browser.