Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 57

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

1file ->several

File size: 13.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 //Trigger information
96 Trigger *TRIG = new Trigger();
97 TRIG->TriggerReader("data/trigger.dat");
98
99 //Propagation of tracks in the B field
100 TrackPropagation *TRACP = new TrackPropagation();
101
102 //Jet information
103 JetsUtil *JETRUN = new JetsUtil();
104
105 //VFD information
106 VeryForward * VFD = new VeryForward();
107
108 //todo(LogName.c_str());
109
110 DataConverter *converter=0;
111
112 if(strstr(line.c_str(),".hep"))
113 {
114 cout<<"#**********************************************************************"<<endl;
115 cout<<"#********** StdHEP file format detected *************"<<endl;
116 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
117 cout<<"#**********************************************************************"<<endl;
118 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
119 }
120 else if(strstr(line.c_str(),".lhe"))
121 {
122 cout<<"#**********************************************************************"<<endl;
123 cout<<"#*********** LHEF file format detected ************"<<endl;
124 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
125 cout<<"#**********************************************************************"<<endl;
126 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
127 }
128 else if(strstr(line.c_str(),".root"))
129 {
130 cout<<"#**********************************************************************"<<endl;
131 cout<<"#********** h2root file format detected *************"<<endl;
132 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
133 cout<<"#**********************************************************************"<<endl;
134 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
135 }
136 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
137
138 TChain chain("GEN");
139 chain.Add(outputfilename.c_str());
140 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
141 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
142 TIter itGen((TCollection*)branchGen);
143
144 //write the output root file
145 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
146 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
147 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
148 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
149 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
150 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
151 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
152 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
153 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
154 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
155 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
156 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
157
158
159 TRootGenParticle *particle;
160 TRootETmis *elementEtmis;
161 TRootElectron *elementElec;
162 TRootMuon *elementMu;
163 TRootPhoton *elementPhoton;
164 TRootTracks *elementTracks;
165 TRootCalo *elementCalo;
166
167 TLorentzVector genMomentum(0,0,0,0);
168 TLorentzVector genMomentumCalo(0,0,0,0);
169 LorentzVector jetMomentum;
170
171 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
172 vector<fastjet::PseudoJet> sorted_jets;
173
174 vector<TLorentzVector> TrackCentral;
175 vector<PhysicsTower> towers;
176
177 vector<TLorentzVector> electron;
178 vector<int> elecPID;
179 vector<TLorentzVector> muon;
180 vector<int> muonPID;
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 elecPID.clear();
198 muonPID.clear();
199 NFCentralQ.Clear();
200
201 itGen.Reset();
202 TrackCentral.clear();
203 towers.clear();
204 input_particles.clear();
205
206 // Loop over all particles in event
207 while( (particle = (TRootGenParticle*) itGen.Next()) )
208 {
209 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
210
211 int pid = abs(particle->PID);
212 float eta = fabs(genMomentum.Eta());
213
214 //subarray of the quarks (i.e. not final particles, i.e status not equal to 1)
215 // in the tracker (i.e. for those we know the tracks)
216 // with enough p_T
217 //// This subarray is needed for the B-jet algorithm
218 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
219 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 ?
220 eta < DET->MAX_TRACKER &&
221 particle->Status != 1 &&
222 particle->PT > DET->PT_QUARKS_MIN ) {
223 NFCentralQ.Add(particle);
224 }
225
226
227 // keeps only final particles, visible by the central detector, including the fiducial volume
228 // the ordering of conditions have been optimised for speed : put first the STATUS condition
229 //
230 //
231 if( (particle->Status == 1) &&
232 (
233 (pid == pMU && eta < DET->MAX_MU) ||
234 (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD)
235 )
236 ) {
237 // TRACP->Propagation(particle,genMomentum);
238 eta = fabs(genMomentum.Eta());
239 switch(pid) {
240
241 case pE: // all electrons with eta < DET->MAX_CALO_FWD
242 DET->SmearElectron(genMomentum);
243 electron.push_back(genMomentum);
244 elecPID.push_back(particle->PID);
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->MAX_TRACKER) {
249 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
250 elementPhoton->Set(genMomentum);
251 }
252 break; // case pGAMMA
253 case pMU: // all muons with eta < DET->MAX_MU
254 DET->SmearMu(genMomentum);
255 muonPID.push_back(particle->PID);
256 muon.push_back(genMomentum);
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
270 if(genMomentum.E() !=0) {
271 if(pid !=pMU) {
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 //genMomentumCalo.SetPtEtaPhiE(CaloTower.Et(),CaloTower.iEta(),CaloTower.iPhi(),CaloTower.fourVector.E);
278 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
279 elementCalo = (TRootCalo*) branchCalo->NewEntry();
280 elementCalo->Set(genMomentumCalo);
281 }
282 }
283
284 // all final charged particles
285 if(
286 ((rand()%100) < DET->TRACKING_EFF) &&
287 (genMomentum.E()!=0) &&
288 (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) &&
289 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
290 (pid != pGAMMA) &&
291 (pid != pPI0) &&
292 (pid != pK0L) &&
293 (pid != pN) &&
294 (pid != pSIGMA0) &&
295 (pid != pDELTA0) &&
296 (pid != pK0S) // not charged particles : invisible by tracker
297 )
298 {
299 elementTracks = (TRootTracks*) branchTracks->NewEntry();
300 elementTracks->Set(genMomentum);
301 TrackCentral.push_back(genMomentum);
302 }
303
304 } // switch
305
306 VFD->ZDC(treeWriter,branchZDC,particle);
307 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
308
309 } // while
310
311 for(unsigned int i=0; i < electron.size(); i++) {
312 if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt)
313 {
314 elementElec = (TRootElectron*) branchElectron->NewEntry();
315 elementElec->Set(electron[i]);
316 elementElec->Charge = sign(elecPID[i]);
317 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
318 }
319 }
320 for(unsigned int i=0; i < muon.size(); i++) {
321 if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt)
322 {
323 elementMu = (TRootMuon*) branchMuon->NewEntry();
324 elementMu->Charge = sign(muonPID[i]);
325 elementMu->Set(muon[i]);
326 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
327 }
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 sorted_jets=JETRUN->RunJets(input_particles);
345 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
346 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
347
348
349 treeWriter->Fill();
350
351 // Add here the trigger
352 // Should test all the trigger table on the event, based on reconstructed objects
353
354 // Trigger *TRIG = new Trigger();
355 // TRIG->TriggerReader("data/trigger.dat");
356
357 } // Loop over all events
358
359 treeWriter->Write();
360
361 cout << "** Exiting..." << endl;
362
363 delete treeWriter;
364 delete treeReader;
365 delete DET;
366 if(converter) delete converter;
367
368 todo("TODO");
369}
370
Note: See TracBrowser for help on using the repository browser.