Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 262

Last change on this file since 262 was 260, checked in by severine ovyn, 15 years ago

add header

File size: 30.3 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \----------------------------------------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
12** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
13** FROG: [hep-ex/0901.2718v1] **
14** **
15** ------------------------------------------------------------------ **
16** **
17** Main authors: **
18** ------------- **
19** **
20** Severine Ovyn Xavier Rouby **
21** severine.ovyn@uclouvain.be xavier.rouby@cern **
22** **
23** Center for Particle Physics and Phenomenology (CP3) **
24** Universite catholique de Louvain (UCL) **
25** Louvain-la-Neuve, Belgium **
26** **
27** Copyright (C) 2008-2009, **
28** All rights reserved. **
29** **
30***********************************************************************/
31
32/// \file Delphes.cpp
33/// \brief executable for the Delphes
34
35#include "TChain.h"
36#include "TApplication.h"
37#include "TStopwatch.h"
38#include "TFile.h"
39
40#include "ExRootTreeReader.h"
41#include "ExRootTreeWriter.h"
42#include "ExRootTreeBranch.h"
43
44#include "DataConverter.h"
45#include "HEPEVTConverter.h"
46#include "LHEFConverter.h"
47#include "STDHEPConverter.h"
48
49#include "SmearUtil.h"
50#include "BFieldProp.h"
51#include "TriggerUtil.h"
52#include "VeryForward.h"
53#include "JetsUtil.h"
54#include "FrogUtil.h"
55
56#include <vector>
57#include <iostream>
58
59#include "Utilities/ExRootAnalysis/interface/ExRootProgressBar.h"
60
61using namespace std;
62
63//------------------------------------------------------------------------------
64void todo(string filename) {
65 ifstream infile(filename.c_str());
66 cout << "** TODO list ..." << endl;
67 while(infile.good()) {
68 string temp;
69 getline(infile,temp);
70 cout << "*" << temp << endl;
71 }
72 cout << "** done...\n";
73}
74
75//------------------------------------------------------------------------------
76
77int main(int argc, char *argv[])
78{
79 int appargc = 2;
80 char *appName= new char[20];
81 char *appOpt= new char[20];
82 sprintf(appName,"Delphes");
83 sprintf(appOpt,"-b");
84 char *appargv[] = {appName,appOpt};
85 TApplication app(appName, &appargc, appargv);
86 delete [] appName;
87 delete [] appOpt;
88
89
90 if(argc != 3 && argc != 4 && argc != 5) {
91 cout << " Usage: " << argv[0] << " input_file output_file [detector_card] [trigger_card] " << endl;
92 cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl;
93 cout << " output_file - output file." << endl;
94 cout << " detector_card - Datacard containing resolution variables for the detector simulation (optional) "<<endl;
95 cout << " trigger_card - Datacard containing the trigger algorithms (optional) "<<endl;
96 exit(1);
97 }
98
99
100
101
102cout << endl << endl;
103
104cout <<"*********************************************************************"<< endl;
105cout <<"*********************************************************************"<< endl;
106cout <<"** **"<< endl;
107cout <<"** Welcome to **"<< endl;
108cout <<"** **"<< endl;
109cout <<"** **"<< endl;
110cout <<"** .ddddddd- lL hH **"<< endl;
111cout <<"** -Dd` `dD: Ll hH` **"<< endl;
112cout <<"** dDd dDd eeee. lL .pp+pp Hh+hhh` -eeee- `sssss **"<< endl;
113cout <<"** -Dd `DD ee. ee Ll .Pp. PP Hh. HH. ee. ee sSs **"<< endl;
114cout <<"** dD` dDd eEeee: lL. pP. pP hH hH` eEeee:` -sSSSs. **"<< endl;
115cout <<"** .Dd :dd eE. LlL PpppPP Hh Hh eE sSS **"<< endl;
116cout <<"** dddddd:. eee+: lL. pp. hh. hh eee+ sssssS **"<< endl;
117cout <<"** Pp **"<< endl;
118cout <<"** **"<< endl;
119cout <<"** Delphes, a framework for the fast simulation **"<< endl;
120cout <<"** of a generic collider experiment **"<< endl;
121cout <<"** **"<< endl;
122cout <<"** --- Version 1.4beta of Delphes --- **"<< endl;
123cout <<"** Last date of change: 9 February 2009 **"<< endl;
124cout <<"** **"<< endl;
125cout <<"** **"<< endl;
126cout <<"** This package uses: **"<< endl;
127cout <<"** ------------------ **"<< endl;
128cout <<"** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **"<< endl;
129cout <<"** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **"<< endl;
130cout <<"** FROG: [hep-ex/0901.2718v1] **"<< endl;
131cout <<"** **"<< endl;
132cout <<"**-----------------------------------------------------------------**"<< endl;
133cout <<"** **"<< endl;
134cout <<"** Main authors: **"<< endl;
135cout <<"** ------------- **"<< endl;
136cout <<"** **"<< endl;
137cout <<"** Séverine Ovyn Xavier Rouby **"<< endl;
138cout <<"** severine.ovyn@uclouvain.be xavier.rouby@cern **"<< endl;
139cout <<"** Center for Particle Physics and Phenomenology (CP3) **"<< endl;
140cout <<"** Universite Catholique de Louvain (UCL) **"<< endl;
141cout <<"** Louvain-la-Neuve, Belgium **"<< endl;
142cout <<"** **"<< endl;
143cout <<"**-----------------------------------------------------------------**"<< endl;
144cout <<"** **"<< endl;
145cout <<"** Former Delphes versions and documentation can be found on : **"<< endl;
146cout <<"** http://www.fynu.ucl.ac.be/delphes.html **"<< endl;
147cout <<"** **"<< endl;
148cout <<"** **"<< endl;
149cout <<"** Disclaimer: this program is a beta version of Delphes and **"<< endl;
150cout <<"** therefore comes without guarantees. Beware of errors and please **"<< endl;
151cout <<"** give us your feedbacks about potential bugs **"<< endl;
152cout <<"** **"<< endl;
153cout <<"*********************************************************************"<< endl;
154cout <<"*********************************************************************"<< endl;
155
156// 1. ********** initialisation ***********
157
158 srand (time (NULL)); /* Initialisation du générateur */
159 TStopwatch globalwatch, loopwatch, triggerwatch, frogwatch;
160 globalwatch.Start();
161
162
163 //read the output TROOT file
164 string inputFileList(argv[1]), outputfilename(argv[2]);
165 if(outputfilename.find(".root") > outputfilename.length()) {
166 cout <<"** ERROR: 'output_file' should be a .root file. Exiting... **"<< endl;
167 exit(1);
168 }
169 //create output log-file name
170 string forLog = outputfilename;
171 string LogName = forLog.erase(forLog.find(".root"));
172 LogName = LogName+"_run.log";
173
174 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
175 outputFile->Close();
176
177 string line;
178 ifstream infile(inputFileList.c_str());
179 infile >> line; // the first line determines the type of input files
180
181 //read the datacard input file
182 string DetDatacard(""); //for detector smearing parameters
183 string TrigDatacard("data/trigger.dat"); //for trigger selection
184
185 string lineCard1,lineCard2;
186 bool detecCard=false,trigCard=false;
187 if(argv[3])
188 {
189 ifstream infile1(argv[3]);
190 infile1 >> lineCard1; // the first line determines the type of input files
191 if(strstr(lineCard1.c_str(),"DETECTOR") && detecCard==true)
192 cerr <<"** ERROR: A DETECTOR card has already been loaded **"<< endl;
193 else if(strstr(lineCard1.c_str(),"DETECTOR") && detecCard==false){DetDatacard =argv[3]; detecCard=true;}
194 else if(strstr(lineCard1.c_str(),"TRIGGER") && trigCard==true)
195 cerr <<"** ERROR: A TRIGGER card has already been loaded **"<< endl;
196 else if(strstr(lineCard1.c_str(),"TRIGGER") && trigCard==false){TrigDatacard =argv[3]; trigCard=true;}
197 }
198 if(argv[4])
199 {
200 ifstream infile2(argv[4]);
201 infile2 >> lineCard2; // the first line determines the type of input files
202 if(strstr(lineCard2.c_str(),"DETECTOR") && detecCard==true)
203 cerr <<"** ERROR: A DETECTOR card has already been loaded **"<< endl;
204 else if(strstr(lineCard2.c_str(),"DETECTOR") && detecCard==false){DetDatacard =argv[4]; detecCard=true;}
205 else if(strstr(lineCard2.c_str(),"TRIGGER") && trigCard==true)
206 cerr <<"** ERROR: A TRIGGER card has already been loaded **"<< endl;
207 else if(strstr(lineCard2.c_str(),"TRIGGER") && trigCard==false){TrigDatacard =argv[4]; trigCard=true;}
208 }
209
210 //Smearing information
211 RESOLution *DET = new RESOLution();
212 cout <<"** **"<< endl;
213 cout <<"** ####### Start reading DETECTOR parameters ####### **"<< endl;
214 cout << left << setw(40) <<"** Opening configuration card: "<<""
215 << left << setw(27) << DetDatacard <<""
216 << right << setw(2) <<"**"<<""<<endl;
217 DET->ReadDataCard(DetDatacard);
218 DET->Logfile(LogName);
219 cout << left << setw(42) <<"** Parameters summarised in file: "<<""
220 << left << setw(27) << LogName <<""
221 << right << setw(2) <<"**"<<""<<endl;
222 cout <<"** **"<< endl;
223
224 //Trigger information
225 cout <<"** ########### Start reading TRIGGER card ########## **"<< endl;
226 if(trigCard==false)
227 {
228 cout <<"** WARNING: Datadard not found, use default card **" << endl;
229 TrigDatacard="data/trigger.dat";
230 }
231 TriggerTable *TRIGT = new TriggerTable();
232 TRIGT->TriggerCardReader(TrigDatacard.c_str());
233 TRIGT->PrintTriggerTable(LogName);
234 if(DET->FLAG_trigger == 1)
235 {
236 cout << left << setw(40) <<"** Opening configuration card: "<<""
237 << left << setw(27) << TrigDatacard <<""
238 << right << setw(2) <<"**"<<""<<endl;
239 cout <<"** **"<< endl;
240 }
241
242 //Propagation of tracks in the B field
243 TrackPropagation *TRACP = new TrackPropagation(DET);
244 //TrackPropagation *TRACP = new TrackPropagation(DetDatacard);
245
246 //Jet information
247 JetsUtil *JETRUN = new JetsUtil(DET);
248 //JetsUtil *JETRUN = new JetsUtil(DetDatacard);
249
250 //VFD information
251 VeryForward * VFD = new VeryForward(DET);
252 //VeryForward * VFD = new VeryForward(DetDatacard);
253
254 // data converters
255 DataConverter *converter=NULL;
256 cout <<"** **"<<endl;
257 cout <<"** ####### Start convertion to TRoot format ######## **"<< endl;
258
259 if(strstr(line.c_str(),".hep"))
260 {
261 cout <<"** StdHEP file format detected **"<<endl;
262 cout <<"** This can take several minutes **"<< endl;
263 converter = new STDHEPConverter(inputFileList,outputfilename);//case ntpl file in input list
264 }
265 else if(strstr(line.c_str(),".lhe"))
266 {
267 cout <<"** LHEF file format detected **"<<endl;
268 cout <<"** This can take several minutes **"<< endl;
269 converter = new LHEFConverter(inputFileList,outputfilename);//case ntpl file in input list
270 }
271 else if(strstr(line.c_str(),".root"))
272 {
273 cout <<"** h2root file format detected **"<<endl;
274 cout <<"** This can take several minutes **"<< endl;
275 converter = new HEPEVTConverter(inputFileList,outputfilename);//case ntpl file in input list
276 }
277 else {
278 cerr << left << setw(4) <<"** "<<""
279 << left << setw(63) << line.c_str() <<""
280 << right << setw(2) <<"**"<<endl;
281 cerr <<"** ERROR: File format not identified -- Exiting... **"<< endl;
282 cout <<"** **"<< endl;
283 cout <<"*********************************************************************"<< endl;
284 return -1;};
285 cout <<"** Exiting conversion... **"<< endl;
286
287 TChain chain("GEN");
288 chain.Add(outputfilename.c_str());
289 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
290 const TClonesArray *branchGen = treeReader->UseBranch("Particle");
291
292 TIter itGen((TCollection*)branchGen);
293
294 //Output file : contents of the analysis object data
295 ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis");
296 ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class());
297 ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class());
298 ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class());
299 ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class());
300 ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class());
301 ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class());
302 ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class());
303 ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class());
304 ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class());
305 ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class());
306 ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class());
307
308 TRootGenParticle *particle;
309 TRootETmis *elementEtmis;
310 TRootElectron *elementElec;
311 TRootMuon *elementMu;
312 TRootPhoton *elementPhoton;
313 TRootTracks *elementTracks;
314 TRootCalo *elementCalo;
315
316 TLorentzVector genMomentum(0,0,0,0); // four-momentum at the vertex
317 TLorentzVector genMomentumBfield(0,0,0,0); // four-momentum at the exit of the tracks
318 TLorentzVector momentumCaloSegmentation(0,0,0,0); // four-momentum in the calo, after applying the calo segmentation
319 LorentzVector jetMomentum;
320
321 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
322 vector<fastjet::PseudoJet> sorted_jets;
323 vector<TLorentzVector> TrackCentral;
324 vector<PhysicsTower> towers;
325 vector<ParticleUtil> electron;
326 vector<ParticleUtil> muon;
327 vector<ParticleUtil> gamma;
328
329 TSimpleArray<TRootGenParticle> NFCentralQ;
330 float iPhi=0,iEta=0;
331
332
333
334// 2. ********** Loop over all events ***********
335 Long64_t entry, allEntries = treeReader->GetEntries();
336 cout <<"** **"<<endl;
337 cout <<"** ####### Start fast detector simulation ######## **"<< endl;
338 cout << left << setw(52) <<"** Total number of events to run: "<<""
339 << left << setw(15) << allEntries <<""
340 << right << setw(2) <<"**"<<endl;
341
342 ExRootProgressBar *Progress = new ExRootProgressBar(allEntries);
343
344 loopwatch.Start();
345
346 // loop on all events
347 for(entry = 0; entry < allEntries; ++entry)
348 {
349 Progress->Update(entry);
350 TLorentzVector PTmis(0,0,0,0);
351 treeReader->ReadEntry(entry);
352 treeWriter->Clear();
353
354 electron.clear();
355 muon.clear();
356 gamma.clear();
357 NFCentralQ.Clear();
358
359 TrackCentral.clear();
360 towers.clear();
361 input_particles.clear();
362
363 // 2.1 Loop over all particles in event
364 itGen.Reset();
365 while( (particle = (TRootGenParticle*) itGen.Next()) ) {
366 int pid = abs(particle->PID);
367
368
369 // 2.1.1********************* preparation for the b-tagging
370 //// This subarray is needed for the B-jet algorithm
371 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
372 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 ?
373 fabs(particle->Eta) < DET->CEN_max_tracker &&
374 particle->Status != 1 &&
375 particle->PT > DET->PT_QUARKS_MIN ) {
376 NFCentralQ.Add(particle);
377 }
378
379 // 2.1.2 ********************* central detector: keep only visible particles
380 // keeps only final particles, visible by the central detector, including the fiducial volume
381 // the ordering of conditions have been optimised for speed : put first the STATUS condition
382 if( (particle->Status == 1) &&
383 ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) &&
384 (fabs(particle->Eta) < DET->CEN_max_calo_fwd)
385 )
386 {
387 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
388 genMomentumBfield = genMomentum;
389
390 // 2.1.2.1 ********************* central detector: magnetic field
391 // genMomentumBfield is then changed with respect to the magnetic field
392 if(DET->FLAG_bfield==1) TRACP->Propagation(particle,genMomentumBfield);
393 float eta=fabs(genMomentumBfield.Eta());
394
395
396 // 2.1.2.2 ********************* central detector: smearing (calorimeters, muon chambers)
397 switch(pid) {
398
399 case pE: // all electrons with eta < DET->MAX_CALO_FWD
400 DET->SmearElectron(genMomentumBfield);
401 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_elec){
402 electron.push_back(ParticleUtil(genMomentumBfield,particle->PID));
403 }
404 break; // case pE
405 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
406 DET->SmearElectron(genMomentumBfield);
407 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_tracker && genMomentumBfield.Pt() > DET->PTCUT_gamma) {
408 gamma.push_back(ParticleUtil(genMomentumBfield,particle->PID));
409 }
410 break; // case pGAMMA
411 case pMU: // all muons with eta < DET->MAX_MU
412 DET->SmearMu(genMomentumBfield);
413 if(genMomentumBfield.E()!=0 && eta < DET->CEN_max_mu && genMomentumBfield.Pt() > DET->PTCUT_muon){
414 muon.push_back(ParticleUtil(genMomentumBfield,particle->PID));
415 }
416 break; // case pMU
417 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
418 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
419 DET->SmearHadron(genMomentumBfield, 0.7);
420 break; // case hadron
421 default: // all other final particles with eta < DET->MAX_CALO_FWD
422 DET->SmearHadron(genMomentumBfield, 1.0);
423 break;
424 } // 2.1.2.2 switch (pid)
425
426
427 // 2.1.2.3 ********************* central detector: calotowers
428 // all final particles but muons and neutrinos
429 // for calorimetric towers and missing PT
430 int charge=Charge(pid);
431 if(genMomentumBfield.E() !=0 && pid != pMU) {
432 // in case the Bfield is not simulated, checks that charged particles have enough pt to reach the calos
433 if ( !DET->FLAG_bfield && charge!=0 && genMomentumBfield.Pt() <= DET->TRACK_ptmin ) { /* particules do not reach calos */ }
434 else { // particles reach calos
435 // applies the calo segmentation and returns iEta & iPhi
436 DET->BinEtaPhi(genMomentumBfield.Phi(), genMomentumBfield.Eta(), iPhi, iEta);
437 if(iEta != -100 && iPhi != -100) {
438 momentumCaloSegmentation.SetPtEtaPhiE(genMomentumBfield.Pt(),iEta,iPhi,genMomentumBfield.E());
439 elementCalo = (TRootCalo*) branchCalo->NewEntry();
440 elementCalo->Set(momentumCaloSegmentation);
441 PhysicsTower Tower(LorentzVector(momentumCaloSegmentation.Px(),momentumCaloSegmentation.Py(),momentumCaloSegmentation.Pz(),momentumCaloSegmentation.E()));
442 towers.push_back(Tower);
443 } // if iEta != -100
444 } // else : when particles reach the calos
445 } // 2.1.2.3 calotowers
446
447
448 // 2.1.2.4 ********************* central detector: tracks
449 // all final charged particles
450 if(
451 (genMomentumBfield.E()!=0) &&
452 (fabs(genMomentumBfield.Eta()) < DET->CEN_max_tracker) &&
453 (DET->FLAG_bfield || ( !DET->FLAG_bfield && genMomentumBfield.Pt() > DET->TRACK_ptmin )) &&
454 // if bfield not simulated, pt should be high enough to be taken into account
455 ((rand()%100) < DET->TRACK_eff) &&
456 (charge!=0)
457 )
458 {
459 elementTracks = (TRootTracks*) branchTracks->NewEntry();
460 elementTracks->Set(genMomentum); // fills px,py,pz,pt,e,eta,phi at vertex
461 elementTracks->Etaout = genMomentumBfield.Eta();
462 elementTracks->Phiout = genMomentumBfield.Phi();
463 // TODO!!! apply a smearing on the position of the origin of the track
464 // elementTracks->SetPosition(particle->X,particle->Y,particle->Z);
465 // uses the output of the bfield computation : Xout=... Yout=... Zout...
466 // TODO!!! elementTrakcs->SetPositionOut(Xout,Yout,Zout);
467 TrackCentral.push_back(genMomentum); // tracks at vertex!
468 } // 2.1.2.4 tracks
469
470 } // 2.1.2 central detector
471
472 // 2.1.3 ********************* forward detectors: zdc
473 if(DET->FLAG_vfd==1) {
474 VFD->ZDC(treeWriter,branchZDC,particle);
475 VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle);
476 }
477
478 } // 2.1 while : loop on all particles of the event.
479
480
481
482 // 2.2 ********** Output preparation & complex objects ***********
483
484 // 2.2.1 ********************* sorting collections by decreasing pt
485 DET->SortedVector(electron);
486 for(unsigned int i=0; i < electron.size(); i++) {
487 elementElec = (TRootElectron*) branchElectron->NewEntry();
488 elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E());
489 elementElec->Charge = sign(electron[i].PID());
490 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);//isolation based on tracks
491 }
492 DET->SortedVector(muon);
493 for(unsigned int i=0; i < muon.size(); i++) {
494 elementMu = (TRootMuon*) branchMuon->NewEntry();
495 elementMu->Charge = sign(muon[i].PID());
496 elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E());
497 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
498 }
499 DET->SortedVector(gamma);
500 for(unsigned int i=0; i < gamma.size(); i++) {
501 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
502 elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E());
503 }
504
505 // 2.2.2 ************* computes the Missing Transverse Momentum
506 TLorentzVector Att(0.,0.,0.,0.);
507 for(unsigned int i=0; i < towers.size(); i++)
508 {
509 Att.SetPxPyPzE(towers[i].fourVector.px, towers[i].fourVector.py, towers[i].fourVector.pz, towers[i].fourVector.E);
510 if(fabs(Att.Eta()) < DET->CEN_max_calo_fwd)
511 {
512 PTmis = PTmis + Att;
513 // create a fastjet::PseudoJet with these components and put it onto
514 // back of the input_particles vector
515 input_particles.push_back(fastjet::PseudoJet(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E));
516 }
517 }
518 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
519 elementEtmis->ET = (PTmis).Pt();
520 elementEtmis->Phi = (-PTmis).Phi();
521 elementEtmis->Px = (-PTmis).Px();
522 elementEtmis->Py = (-PTmis).Py();
523
524 // 2.2.3 ************* B-tag, tau jets
525 sorted_jets=JETRUN->RunJets(input_particles);
526 JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ);
527 JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral);
528
529 treeWriter->Fill();
530 } // 2. Loop over all events ('for' loop)
531
532 cout <<"** **"<< endl;
533 cout <<"** Exiting detector simulation... **"<< endl;
534
535
536 treeWriter->Write();
537 delete treeWriter;
538 loopwatch.Stop();
539
540
541
542 // 3. ********** Trigger & Frog ***********
543 // 3.1 ************ running the trigger in case the FLAG trigger is put to 1 in the datacard
544 triggerwatch.Start();
545 if(DET->FLAG_trigger == 1)
546 {
547 cout <<"** **"<<endl;
548 cout <<"** ########### Start Trigger selection ########### **"<< endl;
549
550 // input
551 TChain chainT("Analysis");
552 chainT.Add(outputfilename.c_str());
553 ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT);
554
555 // output
556 TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron");
557 TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon");
558 TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet");
559 TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet");
560 TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon");
561 TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis");
562
563 ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger");
564 ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class());
565
566
567 Long64_t entryT, allEntriesT = treeReaderT->GetEntries();
568 // loop on all entries
569 for(entryT = 0; entryT < allEntriesT; ++entryT) {
570 treeWriterT->Clear();
571 treeReaderT->ReadEntry(entryT);
572 TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger);
573 treeWriterT->Fill();
574 } // loop on all entries
575 cout <<"** Exiting trigger simulation... **"<< endl;
576
577 treeWriterT->Write();
578 delete treeWriterT;
579 delete treeReaderT;
580 } // trigger
581 triggerwatch.Stop();
582
583
584 // 3.2 ************** FROG display
585 frogwatch.Start();
586 if(DET->FLAG_frog == 1) {
587 cout <<"** **"<<endl;
588 cout <<"** ################## Start FROG ################# **"<< endl;
589
590 FrogDisplay *FROG = new FrogDisplay(DET);
591 FROG->BuildEvents(outputfilename);
592 FROG->BuildGeom();
593 delete FROG;
594 }
595 frogwatch.Stop();
596
597
598
599
600// 4. ********** End & Exit ***********
601
602 globalwatch.Stop();
603 cout <<"** **"<< endl;
604 cout <<"** ################## Time report ################# **"<< endl;
605 cout << left << setw(32) <<"** Time report for "<<""
606 << left << setw(15) << allEntries <<""
607 << right << setw(22) <<"events **"<<endl;
608 cout <<"** **"<< endl;
609 cout << left << setw(10) <<"**"<<""
610 << left << setw(15) <<"Time (s):"<<""
611 << right << setw(15) <<"CPU"<<""
612 << right << setw(15) <<"Real"<<""
613 << right << setw(14) <<"**"<<endl;
614 cout << left << setw(10) <<"**"<<""
615 << left << setw(15) <<" + Global:"<<""
616 << right << setw(15) <<globalwatch.CpuTime()<<""
617 << right << setw(15) <<globalwatch.RealTime()<<""
618 << right << setw(14) <<"**"<<endl;
619 cout << left << setw(10) <<"**"<<""
620 << left << setw(15) <<" + Events:"<<""
621 << right << setw(15) <<loopwatch.CpuTime()<<""
622 << right << setw(15) <<loopwatch.RealTime()<<""
623 << right << setw(14) <<"**"<<endl;
624 if(DET->FLAG_trigger == 1)
625 {
626 cout << left << setw(10) <<"**"<<""
627 << left << setw(15) <<" + Trigger:"<<""
628 << right << setw(15) <<triggerwatch.CpuTime()<<""
629 << right << setw(15) <<triggerwatch.RealTime()<<""
630 << right << setw(14) <<"**"<<endl;
631 }
632 if(DET->FLAG_frog == 1)
633 {
634 cout << left << setw(10) <<"**"<<""
635 << left << setw(15) <<" + Frog:"<<""
636 << right << setw(15) <<frogwatch.CpuTime()<<""
637 << right << setw(15) <<frogwatch.RealTime()<<""
638 << right << setw(14) <<"**"<<endl;
639
640 }
641
642 cout <<"** **"<< endl;
643 cout <<"** Exiting Delphes ... **"<< endl;
644 cout <<"** **"<< endl;
645 cout <<"*********************************************************************"<< endl;
646 cout <<"*********************************************************************"<< endl;
647
648
649 delete treeReader;
650 delete DET;
651// delete TRIGT;
652 delete TRACP;
653 delete JETRUN;
654 delete VFD;
655 delete converter;
656
657// todo("TODO");
658}
Note: See TracBrowser for help on using the repository browser.