Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 179

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

Cleaning, Comments, and removing the pt cut on tracks and calotowers when bfield simulated

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