Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 189

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

tracks with Eta/Phi and Etaout/Phiout. New genMomentumBfield 4-momentum

File size: 18.0 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); // four-momentum at the vertex
175 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
176 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
177 LorentzVector jetMomentum;
178
179 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
180 vector<fastjet::PseudoJet> sorted_jets;
181 vector<TLorentzVector> TrackCentral;
182 vector<PhysicsTower> towers;
183 vector<ParticleUtil> electron;
184 vector<ParticleUtil> muon;
185 vector<ParticleUtil> gamma;
186
187 TSimpleArray<TRootGenParticle> NFCentralQ;
188 float iPhi=0,iEta=0;
189
190
191
192// 2. ********** Loop over all events ***********
193 Long64_t entry, allEntries = treeReader->GetEntries();
194 cout << "** The input list contains " << allEntries << " events" << endl;
195
196 // loop on all events
197 for(entry = 0; entry < allEntries; ++entry)
198 {
199 TLorentzVector PTmis(0,0,0,0);
200 treeReader->ReadEntry(entry);
201 treeWriter->Clear();
202 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
203
204 electron.clear();
205 muon.clear();
206 gamma.clear();
207 NFCentralQ.Clear();
208
209 TrackCentral.clear();
210 towers.clear();
211 input_particles.clear();
212
213 // 2.1 Loop over all particles in event
214 itGen.Reset();
215 while( (particle = (TRootGenParticle*) itGen.Next()) ) {
216 int pid = abs(particle->PID);
217
218
219 // 2.1.1********************* preparation for the b-tagging
220 //// This subarray is needed for the B-jet algorithm
221 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
222 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 ?
223 fabs(particle->Eta) < DET->CEN_max_tracker &&
224 particle->Status != 1 &&
225 particle->PT > DET->PT_QUARKS_MIN ) {
226 NFCentralQ.Add(particle);
227 }
228
229 // 2.1.2 ********************* central detector: keep only visible particles
230 // keeps only final particles, visible by the central detector, including the fiducial volume
231 // the ordering of conditions have been optimised for speed : put first the STATUS condition
232 if( (particle->Status == 1) &&
233 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
234 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
235 )
236 {
237 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
238 genMomentumBfield = genMomentum;
239
240 // 2.1.2.1 ********************* central detector: magnetic field
241 // genMomentumBfield is then changed with respect to the magnetic field
242 if(DET->FLAG_bfield==1) TRACP->Propagation(particle,genMomentumBfield);
243 float eta=fabs(genMomentumBfield.Eta());
244
245
246 // 2.1.2.2 ********************* central detector: smearing (calorimeters, muon chambers)
247 switch(pid) {
248
249 case pE: // all electrons with eta < DET->MAX_CALO_FWD
250 DET->SmearElectron(genMomentumBfield);
251 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_elec){
252 electron.push_back(ParticleUtil(genMomentumBfield,particle->PID));
253 }
254 break; // case pE
255 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
256 DET->SmearElectron(genMomentumBfield);
257 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_gamma) {
258 gamma.push_back(ParticleUtil(genMomentumBfield,particle->PID));
259 }
260 break; // case pGAMMA
261 case pMU: // all muons with eta < DET->MAX_MU
262 DET->SmearMu(genMomentumBfield);
263 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_mu && genMomentumBfield.Pt() > DET->PTCUT_muon){
264 muon.push_back(ParticleUtil(genMomentumBfield,particle->PID));
265 }
266 break; // case pMU
267 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
268 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
269 DET->SmearHadron(genMomentumBfield, 0.7);
270 break; // case hadron
271 default: // all other final particles with eta < DET->MAX_CALO_FWD
272 DET->SmearHadron(genMomentumBfield, 1.0);
273 break;
274 } // 2.1.2.2 switch (pid)
275
276
277 // 2.1.2.3 ********************* central detector: calotowers
278 // all final particles but muons and neutrinos
279 // for calorimetric towers and missing PT
280 int charge=Charge(pid);
281 if(genMomentumBfield.E() !=0 && pid != pMU) {
282 // in case the Bfield is not simulated, checks that charged particles have enough pt to reach the calos
283 if ( !DET->FLAG_bfield && charge!=0 && genMomentumBfield.Pt() <= DET->TRACK_ptmin ) { /* particules do not reach calos */ }
284 else { // particles reach calos
285 // applies the calo segmentation and returns iEta & iPhi
286 DET->BinEtaPhi(genMomentumBfield.Phi(), genMomentumBfield.Eta(), iPhi, iEta);
287 if(iEta != -100 && iPhi != -100) {
288 momentumCaloSegmentation.SetPtEtaPhiE(genMomentumBfield.Pt(),iEta,iPhi,genMomentumBfield.E());
289 elementCalo = (TRootCalo*) branchCalo->NewEntry();
290 elementCalo->Set(momentumCaloSegmentation);
291 PhysicsTower Tower(LorentzVector(momentumCaloSegmentation.Px(),momentumCaloSegmentation.Py(),momentumCaloSegmentation.Pz(),momentumCaloSegmentation.E()));
292 towers.push_back(Tower);
293 } // if iEta != -100
294 } // else : when particles reach the calos
295 } // 2.1.2.3 calotowers
296
297
298 // 2.1.2.4 ********************* central detector: tracks
299 // all final charged particles
300 if(
301 (genMomentumBfield.E()!=0) &&
302 (fabs(genMomentumBfield.Eta()) < DET->CEN_max_tracker) &&
303 (DET->FLAG_bfield || ( !DET->FLAG_bfield && genMomentumBfield.Pt() > DET->TRACK_ptmin )) &&
304 // if bfield not simulated, pt should be high enough to be taken into account
305 ((rand()%100) < DET->TRACK_eff) &&
306 (charge!=0)
307 )
308 {
309 elementTracks = (TRootTracks*) branchTracks->NewEntry();
310 elementTracks->Set(genMomentum); // fills px,py,pz,pt,e,eta,phi at vertex
311 elementTracks->Etaout = genMomentumBfield.Eta();
312 elementTracks->Phiout = genMomentumBfield.Phi();
313 // TODO!!! apply a smearing on the position of the origin of the track
314 // elementTracks->SetPosition(particle->X,particle->Y,particle->Z);
315 // uses the output of the bfield computation : Xout=... Yout=... Zout...
316 // TODO!!! elementTrakcs->SetPositionOut(Xout,Yout,Zout);
317 TrackCentral.push_back(genMomentum); // tracks at vertex!
318 } // 2.1.2.4 tracks
319
320 } // 2.1.2 central detector
321
322 // 2.1.3 ********************* forward detectors: zdc
323 if(DET->FLAG_vfd==1) {
324 VFD->ZDC(treeWriter,branchZDC,particle);
325 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
326 }
327
328 } // 2.1 while : loop on all particles of the event.
329
330
331
332 // 2.2 ********** Output preparation & complex objects ***********
333
334 // 2.2.1 ********************* sorting collections by decreasing pt
335 DET->SortedVector(electron);
336 for(unsigned int i=0; i < electron.size(); i++) {
337 elementElec = (TRootElectron*) branchElectron->NewEntry();
338 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
339 elementElec->Charge = sign(electron[i].PID());
340 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);//isolation based on tracks
341 }
342 DET->SortedVector(muon);
343 for(unsigned int i=0; i < muon.size(); i++) {
344 elementMu = (TRootMuon*) branchMuon->NewEntry();
345 elementMu->Charge = sign(muon[i].PID());
346 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
347 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
348 }
349 DET->SortedVector(gamma);
350 for(unsigned int i=0; i < gamma.size(); i++) {
351 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
352 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
353 }
354
355 // 2.2.2 ************* computes the Missing Transverse Momentum
356 TLorentzVector Att(0.,0.,0.,0.);
357 for(unsigned int i=0; i < towers.size(); i++)
358 {
359 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
360 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
361 {
362 PTmis = PTmis + Att;
363 // create a fastjet::PseudoJet with these components and put it onto
364 // back of the input_particles vector
365 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
366 }
367 }
368 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
369 elementEtmis->ET = (PTmis).Pt();
370 elementEtmis->Phi = (-PTmis).Phi();
371 elementEtmis->Px = (-PTmis).Px();
372 elementEtmis->Py = (-PTmis).Py();
373
374 // 2.2.3 ************* B-tag, tau jets
375 sorted_jets=JETRUN->RunJets(input_particles);
376 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
377 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
378
379 treeWriter->Fill();
380 } // 2. Loop over all events ('for' loop)
381
382 treeWriter->Write();
383 delete treeWriter;
384
385
386
387// 3. ********** Trigger & Frog ***********
388 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
389 if(DET->FLAG_trigger == 1)
390 {
391 // input
392 TChain chainT("Analysis");
393 chainT.Add(outputfilename.c_str());
394 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
395
396 // output
397 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
398 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
399 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
400 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
401 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
402 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
403
404 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
405 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
406
407
408 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
409 cout << "** Trigger: the 'Analysis' tree contains " << allEntriesT << " events" << endl;
410 // loop on all entries
411 for(entryT = 0; entryT < allEntriesT; ++entryT) {
412 treeWriterT->Clear();
413 treeReaderT->ReadEntry(entryT);
414 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
415 treeWriterT->Fill();
416 } // loop on all entries
417
418 treeWriterT->Write();
419 delete treeWriterT;
420 } // trigger
421
422
423 // 3.2 ************** FROG display
424 if(DET->FLAG_frog == 1)
425 {
426 FrogDisplay *FROG = new FrogDisplay();
427 FROG->BuidEvents(outputfilename,DET->NEvents_Frog);
428 FROG->BuildGeom(DetDatacard);
429 }
430
431
432
433
434// 4. ********** End & Exit ***********
435 cout << "** Exiting..." << endl;
436
437 delete treeReader;
438 delete DET;
439 delete TRIGT;
440 delete TRACP;
441 delete JETRUN;
442 delete VFD;
443 if(converter) delete converter;
444
445 todo("TODO");
446}
Note: See TracBrowser for help on using the repository browser.