Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 6

Last change on this file since 6 was 2, checked in by Xavier Rouby, 16 years ago

first commit

File size: 17.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 "Utilities/CDFCones/interface/JetCluAlgorithm.h"
22#include "Utilities/CDFCones/interface/MidPointAlgorithm.h"
23#include "Utilities/CDFCones/interface/PhysicsTower.h"
24#include "Utilities/CDFCones/interface/Cluster.h"
25
26#include "H_BeamParticle.h"
27#include "H_BeamLine.h"
28#include "H_RomanPot.h"
29
30#include "interface/DataConverter.h"
31#include "interface/HEPEVTConverter.h"
32#include "interface/LHEFConverter.h"
33#include "interface/STDHEPConverter.h"
34
35#include "interface/SmearUtil.h"
36
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;
56 char *appName = "JetsSim";
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 }
76
77 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
78 outputFile->Close();
79
80 string line;
81 ifstream infile(inputFileList.c_str());
82 infile >> line; // the first line determines the type of input files
83
84 DataConverter *converter=0;
85
86 if(strstr(line.c_str(),".hep"))
87 {
88 cout<<"*************************************************************************"<<endl;
89 cout<<"************ StdHEP file format detected **************"<<endl;
90 cout<<"************ Starting convertion to TRoot format **************"<<endl;
91 cout<<"*************************************************************************"<<endl;
92 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
93 }
94 else if(strstr(line.c_str(),".lhe"))
95 {
96 cout<<"*************************************************************************"<<endl;
97 cout<<"************ LHEF file format detected **************"<<endl;
98 cout<<"************ Starting convertion to TRoot format **************"<<endl;
99 cout<<"*************************************************************************"<<endl;
100 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
101 }
102 else if(strstr(line.c_str(),".root"))
103 {
104 cout<<"*************************************************************************"<<endl;
105 cout<<"************ h2root file format detected **************"<<endl;
106 cout<<"************ Starting convertion to TRoot format **************"<<endl;
107 cout<<"*************************************************************************"<<endl;
108 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
109 }
110 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
111
112 TChain chain("GEN");
113 chain.Add(outputfilename.c_str());
114 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
115 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
116 TIter itGen((TCollection*)branchGen);
117
118 //write the output root file
119 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
120 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
121 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
122 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
123 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
124 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
125 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
126 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
127 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
128 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
129 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
130 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
131
132
133 TRootGenParticle *particle;
134 TRootETmis *elementEtmis;
135 TRootElectron *elementElec;
136 TRootMuon *elementMu;
137 TRootPhoton *elementPhoton;
138 TRootJet *elementJet;
139 TRootTauJet *elementTauJet;
140 TRootTracks *elementTracks;
141 TRootCalo *elementCalo;
142 TRootZdcHits *elementZdc;
143 TRootRomanPotHits *elementRP220, *elementFP420;
144
145 //read the datacard input file
146 string DetDatacard("");
147 if(argc==4) DetDatacard =argv[3];
148 RESOLution *DET = new RESOLution();
149 DET->ReadDataCard(DetDatacard);
150
151
152 TLorentzVector genMomentum(0,0,0,0);
153 LorentzVector jetMomentum;
154 vector<TLorentzVector> TrackCentral;
155 vector<PhysicsTower> towers;
156 vector<Cluster> jets;
157
158 //Initialisation of Hector
159 extern bool relative_energy;
160 relative_energy = true; // should always be true
161 extern int kickers_on;
162 kickers_on = 1; // should always be 1
163
164 // user should provide : (1) optics file for each beamline, and IPname,
165 // and offset data (s,x) for optical elements
166 H_BeamLine* beamline1 = new H_BeamLine(1,500.);
167 beamline1->fill("data/LHCB1IR5_v6.500.tfs",1,"IP5");
168 beamline1->offsetElements(120,-0.097);
169 H_RomanPot * rp220_1 = new H_RomanPot("rp220_1",220,2000); // RP 220m, 2mm, beam 1
170 H_RomanPot * rp420_1 = new H_RomanPot("rp420_1",420,4000); // RP 420m, 4mm, beam 1
171 beamline1->add(rp220_1);
172 beamline1->add(rp420_1);
173
174 H_BeamLine* beamline2 = new H_BeamLine(1,500.);
175 beamline2->fill("data/LHCB1IR5_v6.500.tfs",-1,"IP5");
176 beamline2->offsetElements(120,+0.097);
177 H_RomanPot * rp220_2 = new H_RomanPot("rp220_2",220,2000);// RP 220m, 2mm, beam 2
178 H_RomanPot * rp420_2 = new H_RomanPot("rp420_2",420,4000);// RP 420m, 4mm, beam 2
179 beamline2->add(rp220_2);
180 beamline2->add(rp420_2);
181
182
183
184 // Loop over all events
185 Long64_t entry, allEntries = treeReader->GetEntries();
186 cout << "** Chain contains " << allEntries << " events" << endl;
187 for(entry = 0; entry < allEntries; ++entry)
188 {
189 TLorentzVector PTmis(0,0,0,0);
190 treeReader->ReadEntry(entry);
191 treeWriter->Clear();
192
193 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
194
195 TSimpleArray<TRootGenParticle> bGen;
196 itGen.Reset();
197 TrackCentral.clear();
198 towers.clear();
199 TSimpleArray<TRootGenParticle> NFCentralQ;
200
201
202 // Loop over all particles in event
203 while( (particle = (TRootGenParticle*) itGen.Next()) )
204 {
205 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
206
207 int pid = abs(particle->PID);
208 float eta = fabs(particle->Eta);
209
210 //subarray of the quarks (i.e. not final particles, i.e status not equal to 1)
211 // in the tracker (i.e. for those we know the tracks)
212 // with enough p_T
213 //// This subarray is needed for the B-jet algorithm
214 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
215 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 ?
216 eta < DET->MAX_TRACKER &&
217 particle->Status != 1 &&
218 particle->PT > DET->PT_QUARKS_MIN ) {
219 NFCentralQ.Add(particle);
220 }
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 if( (particle->Status == 1) &&
226 (
227 (pid == pMU && eta < DET->MAX_MU) ||
228 (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD)
229 )
230 ) {
231 switch(pid) {
232
233 case pE: // all electrons with eta < DET->MAX_CALO_FWD
234 DET->SmearElectron(genMomentum);
235 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
236 elementElec = (TRootElectron*) branchElectron->NewEntry();
237 elementElec->Set(genMomentum);
238 elementElec->Charge = sign(particle->PID);
239 }
240 break; // case pE
241
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
250 case pMU: // all muons with eta < DET->MAX_MU
251 DET->SmearMu(genMomentum);
252 if(genMomentum.E() !=0 ) {
253 elementMu = (TRootMuon*) branchMuon->NewEntry();
254 elementMu->Charge = sign(particle->PID);
255 elementMu->Set(genMomentum);
256 }
257 break; // case pMU
258
259 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
260 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
261 DET->SmearHadron(genMomentum, 0.7);
262 break; // case hadron
263
264 default: // all other final particles with eta < DET->MAX_CALO_FWD
265 DET->SmearHadron(genMomentum, 1.0);
266 break;
267 } // switch (pid)
268
269 // all final particles but muons and neutrinos
270 // for calorimetric towers and mission PT
271 if(genMomentum.E()!=0) {
272 PTmis = PTmis + genMomentum;//ptmis
273 if(pid !=pMU) {
274 towers.push_back(PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E())));
275 elementCalo = (TRootCalo*) branchCalo->NewEntry();
276 elementCalo->Set(genMomentum);
277 }
278 }
279
280
281 // all final charged particles
282 if(
283 ((rand()%100) < DET->TRACKING_EFF) &&
284 (genMomentum.E()!=0) &&
285 (fabs(particle->Eta) < DET->MAX_TRACKER) &&
286 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
287 (pid != pGAMMA) &&
288 (pid != pPI0) &&
289 (pid != pK0L) &&
290 (pid != pN) &&
291 (pid != pSIGMA0) &&
292 (pid != pDELTA0) &&
293 (pid != pK0S) // not charged particles : invisible by tracker
294 )
295 {
296 elementTracks = (TRootTracks*) branchTracks->NewEntry();
297 elementTracks->Set(genMomentum);
298 TrackCentral.push_back(genMomentum);
299 }
300 } // switch
301
302 // Forward particles in CASTOR ?
303/* if (particle->Status == 1 && (fabs(particle->Eta) > DET->MIN_CALO_VFWD)
304 && (fabs(particle->Eta) < DET->MAX_CALO_VFWD)) {
305
306
307 } // CASTOR
308*/
309 // Zero degree calorimeter, for forward neutrons and photons
310 if (particle->Status ==1 && (pid == pN || pid == pGAMMA ) && eta > DET->MIN_ZDC ) {
311// !!!!!!!!! vérifier que particle->Z est bien en micromÚtres!!!
312// !!!!!!!!! vérifier que particle->T est bien en secondes!!!
313// !!!!!!!!! pas de smearing ! on garde trop d'info !
314 elementZdc = (TRootZdcHits*) branchZDC->NewEntry();
315 elementZdc->Set(genMomentum);
316
317 // time of flight t is t = T + d/[ cos(theta) v ]
318 //double tx = acos(particle->Px/particle->Pz);
319 //double ty = acos(particle->Py/particle->Pz);
320 //double theta = (1E-6)*sqrt( pow(tx,2) + pow(ty,2) );
321 //double flight_distance = (DET->ZDC_S - particle->Z*(1E-6))/cos(theta) ; // assumes that Z is in micrometers
322 double flight_distance = (DET->ZDC_S - particle->Z*(1E-6));
323 // assumes also that the emission angle is so small that 1/(cos theta) = 1
324 elementZdc->T = 0*particle->T + flight_distance/speed_of_light; // assumes highly relativistic particles
325 elementZdc->side = sign(particle->Eta);
326
327 }
328
329 // if forward proton
330 if( (pid == pP) && (particle->Status == 1) && (fabs(particle->Eta) > DET->MAX_CALO_FWD) )
331 {
332// !!!!!!!! put here particle->CHARGE and particle->MASS
333 H_BeamParticle p1; /// put here particle->CHARGE and particle->MASS
334 p1.smearAng();
335 p1.smearPos();
336 p1.setPosition(p1.getX()-500.,p1.getY(),p1.getTX()-1*kickers_on*CRANG,p1.getTY(),0);
337 p1.set4Momentum(particle->Px,particle->Py,particle->Pz,particle->E);
338
339 H_BeamLine *beamline;
340 if(particle->Eta >0) beamline = beamline1;
341 else beamline = beamline2;
342
343 p1.computePath(beamline,1);
344
345 if(p1.stopped(beamline)) {
346 if (p1.getStoppingElement()->getName()=="rp220_1" || p1.getStoppingElement()->getName()=="rp220_2") {
347 p1.propagate(DET->RP220_S);
348 elementRP220 = (TRootRomanPotHits*) branchRP220->NewEntry();
349 elementRP220->X = (1E-6)*p1.getX(); // [m]
350 elementRP220->Y = (1E-6)*p1.getY(); // [m]
351 elementRP220->Tx = (1E-6)*p1.getTX(); // [rad]
352 elementRP220->Ty = (1E-6)*p1.getTY(); // [rad]
353 elementRP220->S = p1.getS(); // [m]
354 elementRP220->T = -1; // not yet implemented
355 elementRP220->E = p1.getE(); // not yet implemented
356 elementRP220->q2 = -1; // not yet implemented
357 elementRP220->side = sign(particle->Eta);
358
359 } else if (p1.getStoppingElement()->getName()=="rp420_1" || p1.getStoppingElement()->getName()=="rp420_2") {
360 p1.propagate(DET->FP420_S);
361 elementFP420 = (TRootRomanPotHits*) branchFP420->NewEntry();
362 elementFP420->X = (1E-6)*p1.getX(); // [m]
363 elementFP420->Y = (1E-6)*p1.getY(); // [m]
364 elementFP420->Tx = (1E-6)*p1.getTX(); // [rad]
365 elementFP420->Ty = (1E-6)*p1.getTY(); // [rad]
366 elementFP420->S = p1.getS(); // [m]
367 elementFP420->T = -1; // not yet implemented
368 elementFP420->E = p1.getE(); // not yet implemented
369 elementFP420->q2 = -1; // not yet implemented
370 elementFP420->side = sign(particle->Eta);
371 }
372 }
373
374// if(p1.stopped(beamline) && (p1.getStoppingElement()->getS() > 100))
375// cout << "Eloss =" << 7000.-p1.getE() << " ; " << p1.getStoppingElement()->getName() << endl;
376 } // if forward proton
377
378 } // while
379
380 // computes the Missing Transverse Momentum
381 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
382 elementEtmis->ET = (-PTmis).Pt();
383 elementEtmis->Phi = (-PTmis).Phi();
384 elementEtmis->Px = (-PTmis).Px();
385 elementEtmis->Py = (-PTmis).Py();
386
387 //*****************************
388 jets.clear();
389 switch(DET->JETALGO) {
390 default:
391 case 1: {
392 JetCluAlgorithm jetAlgoC(DET->C_SEEDTHRESHOLD,DET->CONERADIUS,DET->C_ADJACENCYCUT,DET->C_MAXITERATIONS,DET->C_IRATCH,DET->C_OVERLAPTHRESHOLD);
393 jetAlgoC.run(towers, jets);
394 }
395 break;
396
397 case 2: {
398 MidPointAlgorithm jetAlgoM(DET->M_SEEDTHRESHOLD,DET->CONERADIUS,DET->M_CONEAREAFRACTION,DET->M_MAXPAIRSIZE,DET->M_MAXITERATIONS,DET->M_OVERLAPTHRESHOLD);
399 jetAlgoM.run(towers, jets);
400 }
401 break;
402/*
403 case 3: {
404 FastJet
405 }
406 break;
407*/
408 } // switch
409
410 // Loop on all jets
411 // Dealing with jets, tau-jets, b-jets
412 vector<Cluster>::iterator itJet;
413 for(itJet = jets.begin(); itJet != jets.end(); ++itJet) {
414 elementJet = (TRootJet*) branchJet->NewEntry();
415 jetMomentum = itJet->fourVector;
416 TLorentzVector JET;
417 JET.SetPxPyPzE(jetMomentum.px,jetMomentum.py,jetMomentum.pz,jetMomentum.E);
418 elementJet->Set(JET);
419
420 // b-jets
421 bool btag=false;
422 if((fabs(JET.Eta()) < DET->MAX_TRACKER && DET->Btaggedjet(JET, NFCentralQ)))btag=true;
423 elementJet->Btag = btag;
424
425 // Tau jet identification : 1! track and electromagnetic collimation
426 if(fabs(JET.Eta()) < (DET->MAX_TRACKER - DET->TAU_CONE_TRACKS)) {
427 double Energie_tau_central = DET->EnergySmallCone(towers,JET.Eta(),JET.Phi());
428 if(
429 ( Energie_tau_central/JET.E() > DET->TAU_EM_COLLIMATION ) &&
430 ( DET->NumTracks(TrackCentral,DET->PT_TRACK_TAU,JET.Eta(),JET.Phi()) == 1 )
431 ) {
432 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
433 elementTauJet->Set(JET);
434 } // if tau jet
435 } // if JET.eta < tracker - tau_cone : Tau jet identification
436 } // for itJet : loop on all jets
437
438 treeWriter->Fill();
439 // Add here the trigger
440 // Should test all the trigger table on the event, based on reconstructed objects
441 } // Loop over all events
442 treeWriter->Write();
443
444 cout << "** Exiting..." << endl;
445
446 delete treeWriter;
447 delete treeReader;
448 delete DET;
449 if(converter) delete converter;
450
451 todo("TODO");
452}
453
Note: See TracBrowser for help on using the repository browser.