Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 73

Last change on this file since 73 was 73, checked in by uid677, 16 years ago

new PartUtil class

File size: 14.7 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 //Trigger information
97 TriggerTable *TRIGT = new TriggerTable();
98 TRIGT->TriggerCardReader("data/trigger.dat");
99 TRIGT->PrintTriggerTable(LogName);
100
101 //Propagation of tracks in the B field
102 TrackPropagation *TRACP = new TrackPropagation();
103
104 //Jet information
105 JetsUtil *JETRUN = new JetsUtil();
106
107 //VFD information
108 VeryForward * VFD = new VeryForward();
109
110 //todo(LogName.c_str());
111
112 DataConverter *converter=0;
113
114 if(strstr(line.c_str(),".hep"))
115 {
116 cout<<"#**********************************************************************"<<endl;
117 cout<<"#********** StdHEP file format detected *************"<<endl;
118 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
119 cout<<"#**********************************************************************"<<endl;
120 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
121 }
122 else if(strstr(line.c_str(),".lhe"))
123 {
124 cout<<"#**********************************************************************"<<endl;
125 cout<<"#*********** LHEF file format detected ************"<<endl;
126 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
127 cout<<"#**********************************************************************"<<endl;
128 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
129 }
130 else if(strstr(line.c_str(),".root"))
131 {
132 cout<<"#**********************************************************************"<<endl;
133 cout<<"#********** h2root file format detected *************"<<endl;
134 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
135 cout<<"#**********************************************************************"<<endl;
136 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
137 }
138 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
139
140 TChain chain("GEN");
141 chain.Add(outputfilename.c_str());
142 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
143 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
144 TIter itGen((TCollection*)branchGen);
145
146 //write the output root file
147 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
148 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
149 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
150 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
151 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
152 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
153 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
154 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
155 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
156 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
157 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
158 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
159
160
161 TRootGenParticle *particle;
162 TRootETmis *elementEtmis;
163 TRootElectron *elementElec;
164 TRootMuon *elementMu;
165 TRootPhoton *elementPhoton;
166 TRootTracks *elementTracks;
167 TRootCalo *elementCalo;
168
169 TLorentzVector genMomentum(0,0,0,0);
170 TLorentzVector genMomentumCalo(0,0,0,0);
171 LorentzVector jetMomentum;
172
173 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
174 vector<fastjet::PseudoJet> sorted_jets;
175
176 vector<TLorentzVector> TrackCentral;
177 vector<PhysicsTower> towers;
178
179 vector<ParticleUtil> electron;
180 vector<ParticleUtil> muon;
181 TSimpleArray<TRootGenParticle> NFCentralQ;
182
183
184
185 // Loop over all events
186 Long64_t entry, allEntries = treeReader->GetEntries();
187 cout << "** Chain contains " << allEntries << " events" << endl;
188 for(entry = 0; entry < allEntries; ++entry)
189 {
190 TLorentzVector PTmis(0,0,0,0);
191 treeReader->ReadEntry(entry);
192 treeWriter->Clear();
193 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
194
195 electron.clear();
196 muon.clear();
197 NFCentralQ.Clear();
198
199 itGen.Reset();
200 TrackCentral.clear();
201 towers.clear();
202 input_particles.clear();
203
204 // Loop over all particles in event
205 while( (particle = (TRootGenParticle*) itGen.Next()) )
206 {
207 int pid = abs(particle->PID);
208 //// This subarray is needed for the B-jet algorithm
209 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
210 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 ?
211 fabs(particle->Eta) < DET->MAX_TRACKER &&
212 particle->Status != 1 &&
213 particle->PT > DET->PT_QUARKS_MIN ) {
214 NFCentralQ.Add(particle);
215 }
216
217 // keeps only final particles, visible by the central detector, including the fiducial volume
218 // the ordering of conditions have been optimised for speed : put first the STATUS condition
219 //
220 //
221 if( (particle->Status == 1) &&
222 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
223 (fabs(particle->Eta) < DET->MAX_CALO_FWD)
224 )
225 {
226 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
227 TRACP->Propagation(particle,genMomentum);
228 float eta=fabs(genMomentum.Eta());
229
230 switch(pid) {
231
232 case pE: // all electrons with eta < DET->MAX_CALO_FWD
233 DET->SmearElectron(genMomentum);
234 electron.push_back(ParticleUtil(genMomentum,particle->PID));
235 break; // case pE
236 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
237 DET->SmearElectron(genMomentum);
238 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
239 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
240 elementPhoton->Set(genMomentum);
241 }
242 break; // case pGAMMA
243 case pMU: // all muons with eta < DET->MAX_MU
244 DET->SmearMu(genMomentum);
245 muon.push_back(ParticleUtil(genMomentum,particle->PID));
246 break; // case pMU
247 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
248 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
249 DET->SmearHadron(genMomentum, 0.7);
250 break; // case hadron
251 default: // all other final particles with eta < DET->MAX_CALO_FWD
252 DET->SmearHadron(genMomentum, 1.0);
253 break;
254 } // switch (pid)
255
256 // all final particles but muons and neutrinos
257 // for calorimetric towers and mission PT
258
259 if(genMomentum.E() !=0) {
260 if(pid !=pMU) {
261 PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
262 towers.push_back(CaloTower);
263 // create a fastjet::PseudoJet with these components and put it onto
264 // back of the input_particles vector
265 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
266
267 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
268 elementCalo = (TRootCalo*) branchCalo->NewEntry();
269 elementCalo->Set(genMomentumCalo);
270 DET->BinEtaPhi(genMomentumCalo.Phi(), genMomentumCalo.Eta(), elementCalo->Phi, elementCalo->Eta);
271 }
272 }
273
274 // all final charged particles
275 if(
276 ((rand()%100) < DET->TRACKING_EFF) &&
277 (genMomentum.E()!=0) &&
278 (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) &&
279 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
280 (pid != pGAMMA) &&
281 (pid != pPI0) &&
282 (pid != pK0L) &&
283 (pid != pN) &&
284 (pid != pSIGMA0) &&
285 (pid != pDELTA0) &&
286 (pid != pK0S) // not charged particles : invisible by tracker
287 )
288 {
289 elementTracks = (TRootTracks*) branchTracks->NewEntry();
290 elementTracks->Set(genMomentum);
291 TrackCentral.push_back(genMomentum);
292 }
293
294 } // switch
295
296 VFD->ZDC(treeWriter,branchZDC,particle);
297 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
298
299 } // while
300
301 for(unsigned int i=0; i < electron.size(); i++) {
302 if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt)
303 {
304 elementElec = (TRootElectron*) branchElectron->NewEntry();
305 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
306 elementElec->Charge = sign(electron[i].PID());
307 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
308 }
309 }
310 for(unsigned int i=0; i < muon.size(); i++) {
311 if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt)
312 {
313 elementMu = (TRootMuon*) branchMuon->NewEntry();
314 elementMu->Charge = sign(muon[i].PID());
315 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
316 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
317 }
318 }
319
320 // computes the Missing Transverse Momentum
321 TLorentzVector Att(0.,0.,0.,0.);
322 for(unsigned int i=0; i < towers.size(); i++)
323 {
324 Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E);
325 PTmis = PTmis + Att;
326 }
327 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
328 elementEtmis->ET = (PTmis).Pt();
329 elementEtmis->Phi = (-PTmis).Phi();
330 elementEtmis->Px = (-PTmis).Px();
331 elementEtmis->Py = (-PTmis).Py();
332
333 //*****************************
334
335 sorted_jets=JETRUN->RunJets(input_particles);
336 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
337 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
338
339 // Add here the trigger
340 // Should test all the trigger table on the event, based on reconstructed objects
341 treeWriter->Fill();
342
343 } // Loop over all events
344
345 treeWriter->Write();
346 delete treeWriter;
347
348 if(DET->DOTRIGGER == 1)
349 {
350 TChain chainT("Analysis");
351 chainT.Add(outputfilename.c_str());
352 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
353
354 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
355 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
356 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
357 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
358 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
359 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
360
361 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
362 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
363
364 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
365 cout << "** Chain contains " << allEntriesT << " events" << endl;
366 for(entryT = 0; entryT < allEntriesT; ++entryT)
367 {
368 treeWriterT->Clear();
369 treeReaderT->ReadEntry(entryT);
370 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
371 treeWriterT->Fill();
372 }
373
374 treeWriterT->Write();
375 delete treeWriterT;
376 }
377
378 cout << "** Exiting..." << endl;
379
380 delete treeReader;
381 delete DET;
382 if(converter) delete converter;
383
384 todo("TODO");
385}
386
Note: See TracBrowser for help on using the repository browser.