Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 82

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

add default trigger card

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