Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 16

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

remove comments

File size: 19.5 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 "H_BeamParticle.h"
22#include "H_BeamLine.h"
23#include "H_RomanPot.h"
24
25#include "interface/DataConverter.h"
26#include "interface/HEPEVTConverter.h"
27#include "interface/LHEFConverter.h"
28#include "interface/STDHEPConverter.h"
29
30#include "interface/SmearUtil.h"
31#include "Utilities/Fastjet/include/fastjet/PseudoJet.hh"
32#include "Utilities/Fastjet/include/fastjet/ClusterSequence.hh"
33
34// get info on how fastjet was configured
35#include "Utilities/Fastjet/include/fastjet/config.h"
36
37// make sure we have what is needed
38#ifdef ENABLE_PLUGIN_SISCONE
39# include "Utilities/Fastjet/plugins/SISCone/SISConePlugin.hh"
40#endif
41#ifdef ENABLE_PLUGIN_CDFCONES
42# include "Utilities/Fastjet/plugins/CDFCones/CDFMidPointPlugin.hh"
43# include "Utilities/Fastjet/plugins/CDFCones/CDFJetCluPlugin.hh"
44#endif
45
46#include<vector>
47#include<iostream>
48
49
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54void todo(string filename) {
55 ifstream infile(filename.c_str());
56 cout << "** TODO list ..." << endl;
57 while(infile.good()) {
58 string temp;
59 getline(infile,temp);
60 cout << "*" << temp << endl;
61 }
62 cout << "** done...\n";
63}
64
65//------------------------------------------------------------------------------
66
67int main(int argc, char *argv[])
68{
69 int appargc = 2;
70 char *appName = "JetsSim";
71 char *appargv[] = {appName, "-b"};
72 TApplication app(appName, &appargc, appargv);
73
74 if(argc != 4 && argc != 3) {
75 cout << " Usage: " << argv[0] << " input_file" << " output_file" << " data_card " << endl;
76 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
77 cout << " output_file - output file." << endl;
78 cout << " data_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
79 exit(1);
80 }
81
82 srand (time (NULL)); /* Initialisation du générateur */
83
84 //read the input TROOT file
85 string inputFileList(argv[1]), outputfilename(argv[2]);
86 if(outputfilename.find(".root") > outputfilename.length() ) {
87 cout << "output_file should be a .root file!\n";
88 exit(1);
89 }
90
91 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
92 outputFile->Close();
93
94 string line;
95 ifstream infile(inputFileList.c_str());
96 infile >> line; // the first line determines the type of input files
97
98 DataConverter *converter=0;
99
100 if(strstr(line.c_str(),".hep"))
101 {
102 cout<<"*************************************************************************"<<endl;
103 cout<<"************ StdHEP file format detected **************"<<endl;
104 cout<<"************ Starting convertion to TRoot format **************"<<endl;
105 cout<<"*************************************************************************"<<endl;
106 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
107 }
108 else if(strstr(line.c_str(),".lhe"))
109 {
110 cout<<"*************************************************************************"<<endl;
111 cout<<"************ LHEF file format detected **************"<<endl;
112 cout<<"************ Starting convertion to TRoot format **************"<<endl;
113 cout<<"*************************************************************************"<<endl;
114 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
115 }
116 else if(strstr(line.c_str(),".root"))
117 {
118 cout<<"*************************************************************************"<<endl;
119 cout<<"************ h2root file format detected **************"<<endl;
120 cout<<"************ Starting convertion to TRoot format **************"<<endl;
121 cout<<"*************************************************************************"<<endl;
122 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
123 }
124 else { cout << "*** " << line.c_str() << "\n*** file format not identified\n*** Exiting\n"; return -1;};
125
126 TChain chain("GEN");
127 chain.Add(outputfilename.c_str());
128 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
129 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
130 TIter itGen((TCollection*)branchGen);
131
132 //write the output root file
133 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
134 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
135 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
136 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
137 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
138 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
139 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
140 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
141 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
142 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
143 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
144 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
145
146
147 TRootGenParticle *particle;
148 TRootETmis *elementEtmis;
149 TRootElectron *elementElec;
150 TRootMuon *elementMu;
151 TRootPhoton *elementPhoton;
152 TRootJet *elementJet;
153 TRootTauJet *elementTauJet;
154 TRootTracks *elementTracks;
155 TRootCalo *elementCalo;
156 TRootZdcHits *elementZdc;
157 TRootRomanPotHits *elementRP220, *elementFP420;
158
159 //read the datacard input file
160 string DetDatacard("");
161 if(argc==4) DetDatacard =argv[3];
162 RESOLution *DET = new RESOLution();
163 DET->ReadDataCard(DetDatacard);
164
165
166 TLorentzVector genMomentum(0,0,0,0);
167 LorentzVector jetMomentum;
168 vector<TLorentzVector> TrackCentral;
169 vector<PhysicsTower> towers;
170
171 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
172 vector<fastjet::PseudoJet> inclusive_jets;
173 vector<fastjet::PseudoJet> sorted_jets;
174
175 //Initialisation of Hector
176 extern bool relative_energy;
177 relative_energy = true; // should always be true
178 extern int kickers_on;
179 kickers_on = 1; // should always be 1
180
181 // user should provide : (1) optics file for each beamline, and IPname,
182 // and offset data (s,x) for optical elements
183 H_BeamLine* beamline1 = new H_BeamLine(1,500.);
184 beamline1->fill("data/LHCB1IR5_v6.500.tfs",1,"IP5");
185 beamline1->offsetElements(120,-0.097);
186 H_RomanPot * rp220_1 = new H_RomanPot("rp220_1",220,2000); // RP 220m, 2mm, beam 1
187 H_RomanPot * rp420_1 = new H_RomanPot("rp420_1",420,4000); // RP 420m, 4mm, beam 1
188 beamline1->add(rp220_1);
189 beamline1->add(rp420_1);
190
191 H_BeamLine* beamline2 = new H_BeamLine(1,500.);
192 beamline2->fill("data/LHCB1IR5_v6.500.tfs",-1,"IP5");
193 beamline2->offsetElements(120,+0.097);
194 H_RomanPot * rp220_2 = new H_RomanPot("rp220_2",220,2000);// RP 220m, 2mm, beam 2
195 H_RomanPot * rp420_2 = new H_RomanPot("rp420_2",420,4000);// RP 420m, 4mm, beam 2
196 beamline2->add(rp220_2);
197 beamline2->add(rp420_2);
198
199
200
201 // Loop over all events
202 Long64_t entry, allEntries = treeReader->GetEntries();
203 cout << "** Chain contains " << allEntries << " events" << endl;
204 for(entry = 0; entry < allEntries; ++entry)
205 {
206 TLorentzVector PTmis(0,0,0,0);
207 treeReader->ReadEntry(entry);
208 treeWriter->Clear();
209
210 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
211
212 TSimpleArray<TRootGenParticle> bGen;
213 itGen.Reset();
214 TrackCentral.clear();
215 towers.clear();
216 input_particles.clear();
217 TSimpleArray<TRootGenParticle> NFCentralQ;
218 inclusive_jets.clear();
219 sorted_jets.clear();
220
221 // Loop over all particles in event
222 while( (particle = (TRootGenParticle*) itGen.Next()) )
223 {
224 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
225
226 int pid = abs(particle->PID);
227 float eta = fabs(particle->Eta);
228
229 //subarray of the quarks (i.e. not final particles, i.e status not equal to 1)
230 // in the tracker (i.e. for those we know the tracks)
231 // with enough p_T
232 //// This subarray is needed for the B-jet algorithm
233 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
234 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 ?
235 eta < DET->MAX_TRACKER &&
236 particle->Status != 1 &&
237 particle->PT > DET->PT_QUARKS_MIN ) {
238 NFCentralQ.Add(particle);
239 }
240
241
242 // keeps only final particles, visible by the central detector, including the fiducial volume
243 // the ordering of conditions have been optimised for speed : put first the STATUS condition
244 if( (particle->Status == 1) &&
245 (
246 (pid == pMU && eta < DET->MAX_MU) ||
247 (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD)
248 )
249 ) {
250 switch(pid) {
251
252 case pE: // all electrons with eta < DET->MAX_CALO_FWD
253 DET->SmearElectron(genMomentum);
254 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
255 elementElec = (TRootElectron*) branchElectron->NewEntry();
256 elementElec->Set(genMomentum);
257 elementElec->Charge = sign(particle->PID);
258 }
259 break; // case pE
260
261 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
262 DET->SmearElectron(genMomentum);
263 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
264 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
265 elementPhoton->Set(genMomentum);
266 }
267 break; // case pGAMMA
268
269 case pMU: // all muons with eta < DET->MAX_MU
270 DET->SmearMu(genMomentum);
271 if(genMomentum.E() !=0 ) {
272 elementMu = (TRootMuon*) branchMuon->NewEntry();
273 elementMu->Charge = sign(particle->PID);
274 elementMu->Set(genMomentum);
275 }
276 break; // case pMU
277
278 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
279 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
280 DET->SmearHadron(genMomentum, 0.7);
281 break; // case hadron
282
283 default: // all other final particles with eta < DET->MAX_CALO_FWD
284 DET->SmearHadron(genMomentum, 1.0);
285 break;
286 } // switch (pid)
287
288 // all final particles but muons and neutrinos
289 // for calorimetric towers and mission PT
290 if(genMomentum.E()!=0) {
291 PTmis = PTmis + genMomentum;//ptmis
292 if(pid !=pMU) {
293 towers.push_back(PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E())));
294 // create a fastjet::PseudoJet with these components and put it onto
295 // back of the input_particles vector
296 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
297 elementCalo = (TRootCalo*) branchCalo->NewEntry();
298 elementCalo->Set(genMomentum);
299 }
300 }
301
302
303 // all final charged particles
304 if(
305 ((rand()%100) < DET->TRACKING_EFF) &&
306 (genMomentum.E()!=0) &&
307 (fabs(particle->Eta) < DET->MAX_TRACKER) &&
308 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
309 (pid != pGAMMA) &&
310 (pid != pPI0) &&
311 (pid != pK0L) &&
312 (pid != pN) &&
313 (pid != pSIGMA0) &&
314 (pid != pDELTA0) &&
315 (pid != pK0S) // not charged particles : invisible by tracker
316 )
317 {
318 elementTracks = (TRootTracks*) branchTracks->NewEntry();
319 elementTracks->Set(genMomentum);
320 TrackCentral.push_back(genMomentum);
321 }
322 } // switch
323
324 // Forward particles in CASTOR ?
325 /* if (particle->Status == 1 && (fabs(particle->Eta) > DET->MIN_CALO_VFWD)
326 && (fabs(particle->Eta) < DET->MAX_CALO_VFWD)) {
327
328
329 } // CASTOR
330 */
331 // Zero degree calorimeter, for forward neutrons and photons
332 if (particle->Status ==1 && (pid == pN || pid == pGAMMA ) && eta > DET->MIN_ZDC ) {
333 // !!!!!!!!! vérifier que particle->Z est bien en micromÚtres!!!
334 // !!!!!!!!! vérifier que particle->T est bien en secondes!!!
335 // !!!!!!!!! pas de smearing ! on garde trop d'info !
336 elementZdc = (TRootZdcHits*) branchZDC->NewEntry();
337 elementZdc->Set(genMomentum);
338
339 // time of flight t is t = T + d/[ cos(theta) v ]
340 //double tx = acos(particle->Px/particle->Pz);
341 //double ty = acos(particle->Py/particle->Pz);
342 //double theta = (1E-6)*sqrt( pow(tx,2) + pow(ty,2) );
343 //double flight_distance = (DET->ZDC_S - particle->Z*(1E-6))/cos(theta) ; // assumes that Z is in micrometers
344 double flight_distance = (DET->ZDC_S - particle->Z*(1E-6));
345 // assumes also that the emission angle is so small that 1/(cos theta) = 1
346 elementZdc->T = 0*particle->T + flight_distance/speed_of_light; // assumes highly relativistic particles
347 elementZdc->side = sign(particle->Eta);
348
349 }
350
351 // if forward proton
352 if( (pid == pP) && (particle->Status == 1) && (fabs(particle->Eta) > DET->MAX_CALO_FWD) )
353 {
354 // !!!!!!!! put here particle->CHARGE and particle->MASS
355 H_BeamParticle p1; /// put here particle->CHARGE and particle->MASS
356 p1.smearAng();
357 p1.smearPos();
358 p1.setPosition(p1.getX()-500.,p1.getY(),p1.getTX()-1*kickers_on*CRANG,p1.getTY(),0);
359 p1.set4Momentum(particle->Px,particle->Py,particle->Pz,particle->E);
360
361 H_BeamLine *beamline;
362 if(particle->Eta >0) beamline = beamline1;
363 else beamline = beamline2;
364
365 p1.computePath(beamline,1);
366
367 if(p1.stopped(beamline)) {
368 if (p1.getStoppingElement()->getName()=="rp220_1" || p1.getStoppingElement()->getName()=="rp220_2") {
369 p1.propagate(DET->RP220_S);
370 elementRP220 = (TRootRomanPotHits*) branchRP220->NewEntry();
371 elementRP220->X = (1E-6)*p1.getX(); // [m]
372 elementRP220->Y = (1E-6)*p1.getY(); // [m]
373 elementRP220->Tx = (1E-6)*p1.getTX(); // [rad]
374 elementRP220->Ty = (1E-6)*p1.getTY(); // [rad]
375 elementRP220->S = p1.getS(); // [m]
376 elementRP220->T = -1; // not yet implemented
377 elementRP220->E = p1.getE(); // not yet implemented
378 elementRP220->q2 = -1; // not yet implemented
379 elementRP220->side = sign(particle->Eta);
380
381 } else if (p1.getStoppingElement()->getName()=="rp420_1" || p1.getStoppingElement()->getName()=="rp420_2") {
382 p1.propagate(DET->FP420_S);
383 elementFP420 = (TRootRomanPotHits*) branchFP420->NewEntry();
384 elementFP420->X = (1E-6)*p1.getX(); // [m]
385 elementFP420->Y = (1E-6)*p1.getY(); // [m]
386 elementFP420->Tx = (1E-6)*p1.getTX(); // [rad]
387 elementFP420->Ty = (1E-6)*p1.getTY(); // [rad]
388 elementFP420->S = p1.getS(); // [m]
389 elementFP420->T = -1; // not yet implemented
390 elementFP420->E = p1.getE(); // not yet implemented
391 elementFP420->q2 = -1; // not yet implemented
392 elementFP420->side = sign(particle->Eta);
393 }
394 }
395
396 // if(p1.stopped(beamline) && (p1.getStoppingElement()->getS() > 100))
397 // cout << "Eloss =" << 7000.-p1.getE() << " ; " << p1.getStoppingElement()->getName() << endl;
398 } // if forward proton
399
400 } // while
401
402 // computes the Missing Transverse Momentum
403 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
404 elementEtmis->ET = (-PTmis).Pt();
405 elementEtmis->Phi = (-PTmis).Phi();
406 elementEtmis->Px = (-PTmis).Px();
407 elementEtmis->Py = (-PTmis).Py();
408
409 //*****************************
410
411 // we will have four jet definitions, and the first three will be
412 // plugins
413 fastjet::JetDefinition jet_def;
414 fastjet::JetDefinition::Plugin * plugins;
415
416 switch(DET->JETALGO) {
417 default:
418 case 1: {
419
420 // set up a CDF midpoint jet definition
421 #ifdef ENABLE_PLUGIN_CDFCONES
422 plugins = new fastjet::CDFJetCluPlugin(DET->C_SEEDTHRESHOLD,DET->CONERADIUS,DET->C_ADJACENCYCUT,DET->C_MAXITERATIONS,DET->C_IRATCH,DET->C_OVERLAPTHRESHOLD);
423 jet_def = fastjet::JetDefinition(plugins);
424 #else
425 plugins = NULL;
426 #endif
427
428 }
429 break;
430
431 case 2: {
432
433 // set up a CDF midpoint jet definition
434 #ifdef ENABLE_PLUGIN_CDFCONES
435 plugins = new fastjet::CDFMidPointPlugin (DET->M_SEEDTHRESHOLD,DET->CONERADIUS,DET->M_CONEAREAFRACTION,DET->M_MAXPAIRSIZE,DET->M_MAXPAIRSIZE,DET->C_OVERLAPTHRESHOLD);
436 jet_def = fastjet::JetDefinition(plugins);
437 #else
438 plugins = NULL;
439 #endif
440 }
441 break;
442
443 case 3: {
444 // set up a siscone jet definition
445 #ifdef ENABLE_PLUGIN_SISCONE
446 int npass = 0; // do infinite number of passes
447 double protojet_ptmin = 0.0; // use all protojets
448 plugins = new fastjet::SISConePlugin (DET->CONERADIUS,DET->C_OVERLAPTHRESHOLD,npass, protojet_ptmin);
449 jet_def = fastjet::JetDefinition(plugins);
450 #else
451 plugins = NULL;
452 #endif
453 }
454 break;
455
456 case 4: {
457 jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, DET->CONERADIUS);
458 //jet_defs[4] = fastjet::JetDefinition(fastjet::cambridge_algorithm,jet_radius);
459 //jet_defs[5] = fastjet::JetDefinition(fastjet::antikt_algorithm,jet_radius);
460 }
461 break;
462 }
463 // run the jet clustering with the above jet definition
464 if(input_particles.size()!=0)
465 {
466 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
467
468
469 // extract the inclusive jets with pt > 5 GeV
470 double ptmin = 5.0;
471 inclusive_jets = clust_seq.inclusive_jets(ptmin);
472
473 // sort jets into increasing pt
474 sorted_jets = sorted_by_pt(inclusive_jets);
475 }
476 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
477 elementJet = (TRootJet*) branchJet->NewEntry();
478 TLorentzVector JET;
479 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
480 elementJet->Set(JET);
481 // b-jets
482 bool btag=false;
483 if((fabs(JET.Eta()) < DET->MAX_TRACKER && DET->Btaggedjet(JET, NFCentralQ)))btag=true;
484 elementJet->Btag = btag;
485
486 // Tau jet identification : 1! track and electromagnetic collimation
487 if(fabs(JET.Eta()) < (DET->MAX_TRACKER - DET->TAU_CONE_TRACKS)) {
488 double Energie_tau_central = DET->EnergySmallCone(towers,JET.Eta(),JET.Phi());
489 if(
490 ( Energie_tau_central/JET.E() > DET->TAU_EM_COLLIMATION ) &&
491 ( DET->NumTracks(TrackCentral,DET->PT_TRACK_TAU,JET.Eta(),JET.Phi()) == 1 )
492 ) {
493 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
494 elementTauJet->Set(JET);
495 } // if tau jet
496 } // if JET.eta < tracker - tau_cone : Tau jet identification
497 } // for itJet : loop on all jets
498
499 treeWriter->Fill();
500 // Add here the trigger
501 // Should test all the trigger table on the event, based on reconstructed objects
502 } // Loop over all events
503 treeWriter->Write();
504
505 cout << "** Exiting..." << endl;
506
507 delete treeWriter;
508 delete treeReader;
509 delete DET;
510 if(converter) delete converter;
511
512 todo("TODO");
513}
514
Note: See TracBrowser for help on using the repository browser.