Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 98

Last change on this file since 98 was 98, checked in by severine ovyn, 16 years ago

ok for frog

File size: 15.4 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 srand (time (NULL)); /* Initialisation du générateur */
69
70 //read the input TROOT file
71 string inputFileList(argv[1]), outputfilename(argv[2]);
72 if(outputfilename.find(".root") > outputfilename.length() ) {
73 cout << "output_file should be a .root file!\n";
74 exit(1);
75 }
76 //create output log-file name
77 string forLog = outputfilename;
78 string LogName = forLog.erase(forLog.find(".root"));
79 LogName = LogName+"_run.log";
80
81 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
82 outputFile->Close();
83
84 string line;
85 ifstream infile(inputFileList.c_str());
86 infile >> line; // the first line determines the type of input files
87
88 //read the datacard input file
89 string DetDatacard("");
90 if(argc>=4) DetDatacard =argv[3];
91
92 //Smearing information
93 RESOLution *DET = new RESOLution();
94 DET->ReadDataCard(DetDatacard);
95 DET->Logfile(LogName);
96
97 //read the trigger input file
98 string TrigDatacard("data/trigger.dat");
99 if(argc==5) TrigDatacard =argv[4];
100
101 //Trigger information
102 TriggerTable *TRIGT = new TriggerTable();
103 TRIGT->TriggerCardReader(TrigDatacard.c_str());
104 TRIGT->PrintTriggerTable(LogName);
105
106 //Propagation of tracks in the B field
107 TrackPropagation *TRACP = new TrackPropagation();
108
109 //Jet information
110 JetsUtil *JETRUN = new JetsUtil();
111
112 //VFD information
113 VeryForward * VFD = new VeryForward();
114
115 //todo(LogName.c_str());
116
117 DataConverter *converter=0;
118
119 if(strstr(line.c_str(),".hep"))
120 {
121 cout<<"#**********************************************************************"<<endl;
122 cout<<"#********** StdHEP file format detected *************"<<endl;
123 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
124 cout<<"#**********************************************************************"<<endl;
125 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
126 }
127 else if(strstr(line.c_str(),".lhe"))
128 {
129 cout<<"#**********************************************************************"<<endl;
130 cout<<"#*********** LHEF file format detected ************"<<endl;
131 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
132 cout<<"#**********************************************************************"<<endl;
133 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
134 }
135 else if(strstr(line.c_str(),".root"))
136 {
137 cout<<"#**********************************************************************"<<endl;
138 cout<<"#********** h2root file format detected *************"<<endl;
139 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
140 cout<<"#**********************************************************************"<<endl;
141 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
142 }
143 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
144
145 TChain chain("GEN");
146 chain.Add(outputfilename.c_str());
147 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
148 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
149 TIter itGen((TCollection*)branchGen);
150
151 //write the output root file
152 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
153 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
154 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
155 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
156 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
157 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
158 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
159 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
160 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
161 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
162 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
163 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
164
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
181 vector<TLorentzVector> TrackCentral;
182 vector<PhysicsTower> towers;
183
184 vector<ParticleUtil> electron;
185 vector<ParticleUtil> muon;
186 vector<ParticleUtil> gamma;
187
188 TSimpleArray<TRootGenParticle> NFCentralQ;
189
190
191
192 // Loop over all events
193 Long64_t entry, allEntries = treeReader->GetEntries();
194 cout << "** Chain contains " << allEntries << " events" << endl;
195 for(entry = 0; entry < allEntries; ++entry)
196 {
197 TLorentzVector PTmis(0,0,0,0);
198 treeReader->ReadEntry(entry);
199 treeWriter->Clear();
200 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
201
202 electron.clear();
203 muon.clear();
204 gamma.clear();
205 NFCentralQ.Clear();
206
207 TrackCentral.clear();
208 towers.clear();
209 input_particles.clear();
210
211 // Loop over all particles in event
212 itGen.Reset();
213 while( (particle = (TRootGenParticle*) itGen.Next()) )
214 {
215 int pid = abs(particle->PID);
216 //// This subarray is needed for the B-jet algorithm
217 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
218 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 ?
219 fabs(particle->Eta) < DET->CEN_max_tracker &&
220 particle->Status != 1 &&
221 particle->PT > DET->PT_QUARKS_MIN ) {
222 NFCentralQ.Add(particle);
223 }
224
225 // keeps only final particles, visible by the central detector, including the fiducial volume
226 // the ordering of conditions have been optimised for speed : put first the STATUS condition
227 //
228 //
229 if( (particle->Status == 1) &&
230 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
231 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
232 )
233 {
234 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
235 if(DET->FLAG_bfield==1)TRACP->Propagation(particle,genMomentum);
236 float eta=fabs(genMomentum.Eta());
237
238 switch(pid) {
239
240 case pE: // all electrons with eta < DET->MAX_CALO_FWD
241 DET->SmearElectron(genMomentum);
242 if(genMomentum.E()!=0 && eta < DET->CEN_max_tracker && genMomentum.Pt() > DET->PTCUT_elec){
243 electron.push_back(ParticleUtil(genMomentum,particle->PID));
244 }
245 break; // case pE
246 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
247 DET->SmearElectron(genMomentum);
248 if(genMomentum.E()!=0 && eta < DET->CEN_max_tracker && genMomentum.Pt() > DET->PTCUT_gamma) {
249 gamma.push_back(ParticleUtil(genMomentum,particle->PID));
250 }
251 break; // case pGAMMA
252 case pMU: // all muons with eta < DET->MAX_MU
253 DET->SmearMu(genMomentum);
254 if(genMomentum.E()!=0 && eta < DET->CEN_max_mu && genMomentum.Pt() > DET->PTCUT_muon){
255 muon.push_back(ParticleUtil(genMomentum,particle->PID));
256 }
257 break; // case pMU
258 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
259 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
260 DET->SmearHadron(genMomentum, 0.7);
261 break; // case hadron
262 default: // all other final particles with eta < DET->MAX_CALO_FWD
263 DET->SmearHadron(genMomentum, 1.0);
264 break;
265 } // switch (pid)
266
267 // all final particles but muons and neutrinos
268 // for calorimetric towers and mission PT
269 int charge=Charge(pid);
270 if(genMomentum.E() !=0 && pid != pMU) {
271 if(charge == 0 || (charge !=0 && genMomentum.Pt() >= DET->TRACK_ptmin)){
272 PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
273 towers.push_back(CaloTower);
274 // create a fastjet::PseudoJet with these components and put it onto
275 // back of the input_particles vector
276 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
277
278 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
279
280 elementCalo = (TRootCalo*) branchCalo->NewEntry();
281 elementCalo->Set(genMomentumCalo);
282 DET->BinEtaPhi(genMomentumCalo.Phi(), genMomentumCalo.Eta(), elementCalo->Phi, elementCalo->Eta);
283 }
284 }
285
286 // all final charged particles
287 if(
288 (genMomentum.E()!=0) &&
289 (fabs(genMomentum.Eta()) < DET->CEN_max_tracker) &&
290 (genMomentum.Pt() > DET->TRACK_ptmin ) && // pt too small to be taken into account
291 ((rand()%100) < DET->TRACK_eff) &&
292 (charge!=0)
293 )
294 {
295 elementTracks = (TRootTracks*) branchTracks->NewEntry();
296 elementTracks->Set(genMomentum);
297 TrackCentral.push_back(genMomentum);
298 }
299
300 } // switch
301
302 if(DET->FLAG_vfd==1)
303 {
304 VFD->ZDC(treeWriter,branchZDC,particle);
305 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
306 }
307
308 } // while
309
310 DET->SortedVector(electron);
311 for(unsigned int i=0; i < electron.size(); i++) {
312 elementElec = (TRootElectron*) branchElectron->NewEntry();
313 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
314 elementElec->Charge = sign(electron[i].PID());
315 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
316 }
317 DET->SortedVector(muon);
318 for(unsigned int i=0; i < muon.size(); i++) {
319 elementMu = (TRootMuon*) branchMuon->NewEntry();
320 elementMu->Charge = sign(muon[i].PID());
321 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
322 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
323 }
324 DET->SortedVector(gamma);
325 for(unsigned int i=0; i < gamma.size(); i++) {
326 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
327 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
328 }
329
330 // computes the Missing Transverse Momentum
331 TLorentzVector Att(0.,0.,0.,0.);
332 for(unsigned int i=0; i < towers.size(); i++)
333 {
334 Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E);
335 PTmis = PTmis + Att;
336 }
337 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
338 elementEtmis->ET = (PTmis).Pt();
339 elementEtmis->Phi = (-PTmis).Phi();
340 elementEtmis->Px = (-PTmis).Px();
341 elementEtmis->Py = (-PTmis).Py();
342
343 //*****************************
344
345 sorted_jets=JETRUN->RunJets(input_particles);
346 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
347 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
348
349 // Add here the trigger
350 // Should test all the trigger table on the event, based on reconstructed objects
351 treeWriter->Fill();
352
353 } // Loop over all events
354
355 treeWriter->Write();
356 delete treeWriter;
357
358 //running the trigger in case the FLAG trigger is put to 1 in the datacard
359
360 if(DET->FLAG_trigger == 1)
361 {
362 TChain chainT("Analysis");
363 chainT.Add(outputfilename.c_str());
364 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
365
366 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
367 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
368 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
369 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
370 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
371 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
372
373 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
374 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
375
376 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
377 cout << "** Chain contains " << allEntriesT << " events" << endl;
378 for(entryT = 0; entryT < allEntriesT; ++entryT)
379 {
380 treeWriterT->Clear();
381 treeReaderT->ReadEntry(entryT);
382 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
383 treeWriterT->Fill();
384 }
385
386 treeWriterT->Write();
387 delete treeWriterT;
388 }
389
390 //FROG display
391 if(DET->FLAG_frog == 1)
392 {
393 FrogDisplay *FROG = new FrogDisplay();
394 FROG->BuidEvents(outputfilename,DET->NEvents_Frog);
395 FROG->BuildGeom(DetDatacard);
396 }
397
398 cout << "** Exiting..." << endl;
399
400 delete treeReader;
401 delete DET;
402 delete TRIGT;
403 delete TRACP;
404 delete JETRUN;
405 delete VFD;
406
407 if(converter) delete converter;
408
409 todo("TODO");
410}
411
Note: See TracBrowser for help on using the repository browser.