Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 70

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

remove bug trigger

File size: 14.6 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) {
59 cout << " Usage: " << argv[0] << " input_file" << " output_file" << " data_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 << " data_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
63 exit(1);
64 }
65
66 srand (time (NULL)); /* Initialisation du générateur */
67
68 //read the input TROOT file
69 string inputFileList(argv[1]), outputfilename(argv[2]);
70 if(outputfilename.find(".root") > outputfilename.length() ) {
71 cout << "output_file should be a .root file!\n";
72 exit(1);
73 }
74 //create output log-file name
75 string forLog = outputfilename;
76 string LogName = forLog.erase(forLog.find(".root"));
77 LogName = LogName+"_run.log";
78
79 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
80 outputFile->Close();
81
82 string line;
83 ifstream infile(inputFileList.c_str());
84 infile >> line; // the first line determines the type of input files
85
86 //read the datacard input file
87 string DetDatacard("");
88 if(argc==4) DetDatacard =argv[3];
89
90 //Smearing information
91 RESOLution *DET = new RESOLution();
92 DET->ReadDataCard(DetDatacard);
93 DET->Logfile(LogName);
94
95
96 //Trigger information
97 TriggerTable *TRIGT = new TriggerTable();
98 TRIGT->TriggerCardReader("data/trigger.dat");
99 TRIGT->PrintTriggerTable(LogName);
100 //outputFile->Close();
101
102
103 //Propagation of tracks in the B field
104 TrackPropagation *TRACP = new TrackPropagation();
105
106 //Jet information
107 JetsUtil *JETRUN = new JetsUtil();
108
109 //VFD information
110 VeryForward * VFD = new VeryForward();
111
112 //todo(LogName.c_str());
113
114 DataConverter *converter=0;
115
116 if(strstr(line.c_str(),".hep"))
117 {
118 cout<<"#**********************************************************************"<<endl;
119 cout<<"#********** StdHEP file format detected *************"<<endl;
120 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
121 cout<<"#**********************************************************************"<<endl;
122 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
123 }
124 else if(strstr(line.c_str(),".lhe"))
125 {
126 cout<<"#**********************************************************************"<<endl;
127 cout<<"#*********** LHEF file format detected ************"<<endl;
128 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
129 cout<<"#**********************************************************************"<<endl;
130 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
131 }
132 else if(strstr(line.c_str(),".root"))
133 {
134 cout<<"#**********************************************************************"<<endl;
135 cout<<"#********** h2root file format detected *************"<<endl;
136 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
137 cout<<"#**********************************************************************"<<endl;
138 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
139 }
140 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
141
142 TChain chain("GEN");
143 chain.Add(outputfilename.c_str());
144 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
145 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
146 TIter itGen((TCollection*)branchGen);
147
148 //write the output root file
149 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
150 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
151 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
152 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
153 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
154 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
155 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
156 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
157 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
158 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
159 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
160 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
161
162
163 TRootGenParticle *particle;
164 TRootETmis *elementEtmis;
165 TRootElectron *elementElec;
166 TRootMuon *elementMu;
167 TRootPhoton *elementPhoton;
168 TRootTracks *elementTracks;
169 TRootCalo *elementCalo;
170
171 TLorentzVector genMomentum(0,0,0,0);
172 TLorentzVector genMomentumCalo(0,0,0,0);
173 LorentzVector jetMomentum;
174
175 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
176 vector<fastjet::PseudoJet> sorted_jets;
177
178 vector<TLorentzVector> TrackCentral;
179 vector<PhysicsTower> towers;
180
181 vector<TLorentzVector> electron;
182 vector<int> elecPID;
183 vector<TLorentzVector> muon;
184 vector<int> muonPID;
185 TSimpleArray<TRootGenParticle> NFCentralQ;
186
187
188
189 // Loop over all events
190 Long64_t entry, allEntries = treeReader->GetEntries();
191 cout << "** Chain contains " << allEntries << " events" << endl;
192 for(entry = 0; entry < allEntries; ++entry)
193 {
194 TLorentzVector PTmis(0,0,0,0);
195 treeReader->ReadEntry(entry);
196 treeWriter->Clear();
197 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
198
199 electron.clear();
200 muon.clear();
201 elecPID.clear();
202 muonPID.clear();
203 NFCentralQ.Clear();
204
205 itGen.Reset();
206 TrackCentral.clear();
207 towers.clear();
208 input_particles.clear();
209
210 // Loop over all particles in event
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 switch(pid) {
236
237 case pE: // all electrons with eta < DET->MAX_CALO_FWD
238 DET->SmearElectron(genMomentum);
239 electron.push_back(genMomentum);
240 elecPID.push_back(particle->PID);
241 break; // case pE
242 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
243 DET->SmearElectron(genMomentum);
244 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
245 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
246 elementPhoton->Set(genMomentum);
247 }
248 break; // case pGAMMA
249 case pMU: // all muons with eta < DET->MAX_MU
250 DET->SmearMu(genMomentum);
251 muonPID.push_back(particle->PID);
252 muon.push_back(genMomentum);
253 break; // case pMU
254 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
255 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
256 DET->SmearHadron(genMomentum, 0.7);
257 break; // case hadron
258 default: // all other final particles with eta < DET->MAX_CALO_FWD
259 DET->SmearHadron(genMomentum, 1.0);
260 break;
261 } // switch (pid)
262
263 // all final particles but muons and neutrinos
264 // for calorimetric towers and mission PT
265
266 if(genMomentum.E() !=0) {
267 if(pid !=pMU) {
268 PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
269 towers.push_back(CaloTower);
270 // create a fastjet::PseudoJet with these components and put it onto
271 // back of the input_particles vector
272 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
273
274 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
275 elementCalo = (TRootCalo*) branchCalo->NewEntry();
276 elementCalo->Set(genMomentumCalo);
277 }
278 }
279
280 // all final charged particles
281 if(
282 ((rand()%100) < DET->TRACKING_EFF) &&
283 (genMomentum.E()!=0) &&
284 (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) &&
285 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
286 (pid != pGAMMA) &&
287 (pid != pPI0) &&
288 (pid != pK0L) &&
289 (pid != pN) &&
290 (pid != pSIGMA0) &&
291 (pid != pDELTA0) &&
292 (pid != pK0S) // not charged particles : invisible by tracker
293 )
294 {
295 elementTracks = (TRootTracks*) branchTracks->NewEntry();
296 elementTracks->Set(genMomentum);
297 TrackCentral.push_back(genMomentum);
298 }
299
300 } // switch
301
302 VFD->ZDC(treeWriter,branchZDC,particle);
303 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
304
305 } // while
306
307 // computes the Missing Transverse Momentum
308 TLorentzVector Att(0.,0.,0.,0.);
309 for(unsigned int i=0; i < towers.size(); i++)
310 {
311 Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E);
312 PTmis = PTmis + Att;
313 }
314 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
315 elementEtmis->ET = (PTmis).Pt();
316 elementEtmis->Phi = (-PTmis).Phi();
317 elementEtmis->Px = (-PTmis).Px();
318 elementEtmis->Py = (-PTmis).Py();
319 //*****************************
320
321 for(unsigned int i=0; i < electron.size(); i++) {
322 if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt)
323 {
324 elementElec = (TRootElectron*) branchElectron->NewEntry();
325 elementElec->Set(electron[i]);
326 elementElec->Charge = sign(elecPID[i]);
327 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
328 }
329 }
330 for(unsigned int i=0; i < muon.size(); i++) {
331 if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt)
332 {
333 elementMu = (TRootMuon*) branchMuon->NewEntry();
334 elementMu->Charge = sign(muonPID[i]);
335 elementMu->Set(muon[i]);
336 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
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 if(DET->DOTRIGGER == 1)
354 {
355 TChain chainT("Analysis");
356 chainT.Add(outputfilename.c_str());
357 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
358
359 const TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
360 const TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
361 const TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
362 const TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
363 const TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
364 const TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
365
366 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
367 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
368
369 TRootTrigger *elementTrigger;
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(entry);
377 elementTrigger = (TRootTrigger*) branchTrigger->NewEntry();
378 elementTrigger->Accepted=1;
379 treeWriterT->Fill();
380 }
381
382 treeWriterT->Write();
383 delete treeWriterT;
384 }
385
386 cout << "** Exiting..." << endl;
387
388 delete treeReader;
389 delete DET;
390 if(converter) delete converter;
391
392 todo("TODO");
393}
394
Note: See TracBrowser for help on using the repository browser.