Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 63

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

Trigger log file

File size: 13.6 KB
RevLine 
[2]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"
[55]27#include "interface/BFieldProp.h"
28#include "interface/TriggerUtil.h"
29#include "interface/VeryForward.h"
30#include "interface/JetUtils.h"
[2]31
[63]32#include "Utilities/FROG/Examples/Sim_Delphes/FrogUtil.h"
33
[55]34#include <vector>
35#include <iostream>
[11]36
[2]37using namespace std;
38
39//------------------------------------------------------------------------------
40void todo(string filename) {
41 ifstream infile(filename.c_str());
42 cout << "** TODO list ..." << endl;
43 while(infile.good()) {
44 string temp;
45 getline(infile,temp);
46 cout << "*" << temp << endl;
47 }
48 cout << "** done...\n";
49}
50
51//------------------------------------------------------------------------------
52
53int main(int argc, char *argv[])
54{
55 int appargc = 2;
[55]56 char *appName = "Delphes";
[2]57 char *appargv[] = {appName, "-b"};
58 TApplication app(appName, &appargc, appargv);
59
60 if(argc != 4 && argc != 3) {
61 cout << " Usage: " << argv[0] << " input_file" << " output_file" << " data_card " << endl;
62 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
63 cout << " output_file - output file." << endl;
64 cout << " data_card - Datacard containing resolution variables for the detector simulation (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 }
[44]76 //create output log-file name
[45]77 string forLog = outputfilename;
78 string LogName = forLog.erase(forLog.find(".root"));
[44]79 LogName = LogName+"_run.log";
[2]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
[44]87
88 //read the datacard input file
89 string DetDatacard("");
90 if(argc==4) DetDatacard =argv[3];
[55]91
92 //Smearing information
[44]93 RESOLution *DET = new RESOLution();
94 DET->ReadDataCard(DetDatacard);
95 DET->Logfile(LogName);
96
[63]97
[55]98 //Trigger information
[63]99 TriggerTable *TRIGT = new TriggerTable();
100 TRIGT->TriggerCardReader("data/trigger.dat");
101 TRIGT->PrintTriggerTable(LogName);
102 //outputFile->Close();
103
[55]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());
[2]115
116 DataConverter *converter=0;
117
118 if(strstr(line.c_str(),".hep"))
119 {
[44]120 cout<<"#**********************************************************************"<<endl;
121 cout<<"#********** StdHEP file format detected *************"<<endl;
122 cout<<"#*********** Starting convertion to TRoot format **************"<<endl;
123 cout<<"#**********************************************************************"<<endl;
[2]124 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
125 }
126 else if(strstr(line.c_str(),".lhe"))
127 {
[44]128 cout<<"#**********************************************************************"<<endl;
129 cout<<"#*********** LHEF file format detected ************"<<endl;
130 cout<<"#*********** Starting convertion to TRoot format ************"<<endl;
131 cout<<"#**********************************************************************"<<endl;
[2]132 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
133 }
134 else if(strstr(line.c_str(),".root"))
135 {
[44]136 cout<<"#**********************************************************************"<<endl;
137 cout<<"#********** h2root file format detected *************"<<endl;
138 cout<<"#********** Starting convertion to TRoot format *************"<<endl;
139 cout<<"#**********************************************************************"<<endl;
[2]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;};
[30]143
[2]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());
[30]163
164
[2]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);
[38]174 TLorentzVector genMomentumCalo(0,0,0,0);
[2]175 LorentzVector jetMomentum;
[55]176
177 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
178 vector<fastjet::PseudoJet> sorted_jets;
179
[2]180 vector<TLorentzVector> TrackCentral;
181 vector<PhysicsTower> towers;
[11]182
[30]183 vector<TLorentzVector> electron;
184 vector<int> elecPID;
185 vector<TLorentzVector> muon;
186 vector<int> muonPID;
187 TSimpleArray<TRootGenParticle> NFCentralQ;
188
[55]189
190
[30]191 // Loop over all events
[2]192 Long64_t entry, allEntries = treeReader->GetEntries();
193 cout << "** Chain contains " << allEntries << " events" << endl;
194 for(entry = 0; entry < allEntries; ++entry)
195 {
196 TLorentzVector PTmis(0,0,0,0);
197 treeReader->ReadEntry(entry);
198 treeWriter->Clear();
199 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
200
[30]201 electron.clear();
202 muon.clear();
203 elecPID.clear();
204 muonPID.clear();
205 NFCentralQ.Clear();
206
[2]207 itGen.Reset();
208 TrackCentral.clear();
209 towers.clear();
[11]210 input_particles.clear();
[30]211
[2]212 // Loop over all particles in event
213 while( (particle = (TRootGenParticle*) itGen.Next()) )
214 {
[55]215 int pid = abs(particle->PID);
[59]216 //// This subarray is needed for the B-jet algorithm
[2]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 ?
[59]219 fabs(particle->Eta) < DET->MAX_TRACKER &&
[2]220 particle->Status != 1 &&
221 particle->PT > DET->PT_QUARKS_MIN ) {
222 NFCentralQ.Add(particle);
223 }
[59]224
[2]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
[55]227 //
228 //
[2]229 if( (particle->Status == 1) &&
[59]230 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
231 (fabs(particle->Eta) < DET->MAX_CALO_FWD)
232 )
233 {
234 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
[63]235//cout<<"avant "<<genMomentum.Pt()<<endl;
[59]236 TRACP->Propagation(particle,genMomentum);
[63]237//cout<<"apres "<<genMomentum.Pt()<<endl;
[59]238 float eta=fabs(genMomentum.Eta());
[55]239 switch(pid) {
[30]240
241 case pE: // all electrons with eta < DET->MAX_CALO_FWD
242 DET->SmearElectron(genMomentum);
243 electron.push_back(genMomentum);
[63]244 // cout<<"apres le smearing "<<genMomentum.Pt()<<endl;
245
[30]246 elecPID.push_back(particle->PID);
247 break; // case pE
248 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
249 DET->SmearElectron(genMomentum);
250 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
251 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
252 elementPhoton->Set(genMomentum);
253 }
254 break; // case pGAMMA
255 case pMU: // all muons with eta < DET->MAX_MU
256 DET->SmearMu(genMomentum);
257 muonPID.push_back(particle->PID);
258 muon.push_back(genMomentum);
259 break; // case pMU
260 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
261 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
262 DET->SmearHadron(genMomentum, 0.7);
263 break; // case hadron
264 default: // all other final particles with eta < DET->MAX_CALO_FWD
265 DET->SmearHadron(genMomentum, 1.0);
266 break;
267 } // switch (pid)
[59]268
[30]269 // all final particles but muons and neutrinos
270 // for calorimetric towers and mission PT
[43]271
272 if(genMomentum.E() !=0) {
[30]273 if(pid !=pMU) {
[38]274 PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
275 towers.push_back(CaloTower);
[30]276 // create a fastjet::PseudoJet with these components and put it onto
277 // back of the input_particles vector
278 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
[63]279 // cout<<"on remplie avec "<<pid<<" Pt "<<genMomentum.Pt()<<" Eta "<<genMomentum.Eta()<<" Phi "<<genMomentum.Eta()<<endl;
280
[38]281 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
[30]282 elementCalo = (TRootCalo*) branchCalo->NewEntry();
[38]283 elementCalo->Set(genMomentumCalo);
[30]284 }
285 }
286
287 // all final charged particles
288 if(
289 ((rand()%100) < DET->TRACKING_EFF) &&
290 (genMomentum.E()!=0) &&
[55]291 (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) &&
[30]292 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
293 (pid != pGAMMA) &&
294 (pid != pPI0) &&
295 (pid != pK0L) &&
296 (pid != pN) &&
297 (pid != pSIGMA0) &&
298 (pid != pDELTA0) &&
299 (pid != pK0S) // not charged particles : invisible by tracker
300 )
301 {
302 elementTracks = (TRootTracks*) branchTracks->NewEntry();
303 elementTracks->Set(genMomentum);
304 TrackCentral.push_back(genMomentum);
305 }
[59]306
[49]307 } // switch
[30]308
[55]309 VFD->ZDC(treeWriter,branchZDC,particle);
310 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
[11]311
[2]312 } // while
[63]313// cout<<"*************"<<endl;
[30]314 for(unsigned int i=0; i < electron.size(); i++) {
[33]315 if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt)
[30]316 {
317 elementElec = (TRootElectron*) branchElectron->NewEntry();
318 elementElec->Set(electron[i]);
319 elementElec->Charge = sign(elecPID[i]);
320 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
321 }
322 }
323 for(unsigned int i=0; i < muon.size(); i++) {
[33]324 if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt)
[30]325 {
326 elementMu = (TRootMuon*) branchMuon->NewEntry();
327 elementMu->Charge = sign(muonPID[i]);
328 elementMu->Set(muon[i]);
329 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
330 }
331 }
332
[2]333 // computes the Missing Transverse Momentum
[40]334 TLorentzVector Att(0.,0.,0.,0.);
335 for(unsigned int i=0; i < towers.size(); i++)
336 {
337 Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E);
338 PTmis = PTmis + Att;
339 }
[2]340 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
[40]341 elementEtmis->ET = (PTmis).Pt();
[2]342 elementEtmis->Phi = (-PTmis).Phi();
343 elementEtmis->Px = (-PTmis).Px();
344 elementEtmis->Py = (-PTmis).Py();
[59]345
[11]346 //*****************************
[59]347 treeWriter->Fill();
[30]348
[55]349 sorted_jets=JETRUN->RunJets(input_particles);
350 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
351 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
352
[2]353 // Add here the trigger
354 // Should test all the trigger table on the event, based on reconstructed objects
[55]355
[2]356 } // Loop over all events
[55]357
[2]358 treeWriter->Write();
359
360 cout << "** Exiting..." << endl;
361
362 delete treeWriter;
363 delete treeReader;
364 delete DET;
365 if(converter) delete converter;
366
367 todo("TODO");
368}
369
Note: See TracBrowser for help on using the repository browser.