Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 44

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

add log file

File size: 20.9 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 //create output log-file name
91 string LogName = outputfilename.erase(outputfilename.find(".root"));
92 LogName = LogName+"_run.log";
93
94 TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after
95 outputFile->Close();
96
97 string line;
98 ifstream infile(inputFileList.c_str());
99 infile >> line; // the first line determines the type of input files
100
101 //read the datacard input file
102 string DetDatacard("");
103 if(argc==4) DetDatacard =argv[3];
104 RESOLution *DET = new RESOLution();
105 DET->ReadDataCard(DetDatacard);
106 DET->Logfile(LogName);
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 TRootJet *elementJet;
165 TRootTauJet *elementTauJet;
166 TRootTracks *elementTracks;
167 TRootCalo *elementCalo;
168 TRootZdcHits *elementZdc;
169 TRootRomanPotHits *elementRP220, *elementFP420;
170
171 TLorentzVector genMomentum(0,0,0,0);
172 TLorentzVector genMomentumCalo(0,0,0,0);
173 LorentzVector jetMomentum;
174 vector<TLorentzVector> TrackCentral;
175 vector<PhysicsTower> towers;
176 vector<fastjet::PseudoJet> input_particles;//for FastJet algorithm
177 vector<fastjet::PseudoJet> inclusive_jets;
178 vector<fastjet::PseudoJet> sorted_jets;
179
180 vector<TLorentzVector> electron;
181 vector<int> elecPID;
182 vector<TLorentzVector> muon;
183 vector<int> muonPID;
184 TSimpleArray<TRootGenParticle> NFCentralQ;
185
186 //Initialisation of Hector
187 extern bool relative_energy;
188 relative_energy = true; // should always be true
189 extern int kickers_on;
190 kickers_on = 1; // should always be 1
191
192 // user should provide : (1) optics file for each beamline, and IPname,
193 // and offset data (s,x) for optical elements
194 H_BeamLine* beamline1 = new H_BeamLine(1,500.);
195 beamline1->fill("data/LHCB1IR5_v6.500.tfs",1,"IP5");
196 beamline1->offsetElements(120,-0.097);
197 H_RomanPot * rp220_1 = new H_RomanPot("rp220_1",220,2000); // RP 220m, 2mm, beam 1
198 H_RomanPot * rp420_1 = new H_RomanPot("rp420_1",420,4000); // RP 420m, 4mm, beam 1
199 beamline1->add(rp220_1);
200 beamline1->add(rp420_1);
201
202 H_BeamLine* beamline2 = new H_BeamLine(1,500.);
203 beamline2->fill("data/LHCB1IR5_v6.500.tfs",-1,"IP5");
204 beamline2->offsetElements(120,+0.097);
205 H_RomanPot * rp220_2 = new H_RomanPot("rp220_2",220,2000);// RP 220m, 2mm, beam 2
206 H_RomanPot * rp420_2 = new H_RomanPot("rp420_2",420,4000);// RP 420m, 4mm, beam 2
207 beamline2->add(rp220_2);
208 beamline2->add(rp420_2);
209
210 // we will have four jet definitions, and the first three will be
211 // plugins
212 fastjet::JetDefinition jet_def;
213 fastjet::JetDefinition::Plugin * plugins;
214
215 switch(DET->JETALGO) {
216 default:
217 case 1: {
218
219 // set up a CDF midpoint jet definition
220#ifdef ENABLE_PLUGIN_CDFCONES
221 plugins = new fastjet::CDFJetCluPlugin(DET->SEEDTHRESHOLD,DET->CONERADIUS,DET->C_ADJACENCYCUT,DET->C_MAXITERATIONS,DET->C_IRATCH,DET->OVERLAPTHRESHOLD);
222 jet_def = fastjet::JetDefinition(plugins);
223#else
224 plugins = NULL;
225#endif
226 }
227 break;
228
229 case 2: {
230
231 // set up a CDF midpoint jet definition
232#ifdef ENABLE_PLUGIN_CDFCONES
233 plugins = new fastjet::CDFMidPointPlugin (DET->SEEDTHRESHOLD,DET->CONERADIUS,DET->M_CONEAREAFRACTION,DET->M_MAXPAIRSIZE,DET->M_MAXPAIRSIZE,DET->OVERLAPTHRESHOLD);
234 jet_def = fastjet::JetDefinition(plugins);
235#else
236 plugins = NULL;
237#endif
238 }
239 break;
240 case 3: {
241 // set up a siscone jet definition
242#ifdef ENABLE_PLUGIN_SISCONE
243 plugins = new fastjet::SISConePlugin (DET->CONERADIUS,DET->OVERLAPTHRESHOLD,DET->NPASS, DET->PROTOJET_PTMIN);
244 jet_def = fastjet::JetDefinition(plugins);
245#else
246 plugins = NULL;
247#endif
248 }
249 break;
250
251 case 4: {
252 jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, DET->CONERADIUS);
253 //jet_defs[4] = fastjet::JetDefinition(fastjet::cambridge_algorithm,jet_radius);
254 //jet_defs[5] = fastjet::JetDefinition(fastjet::antikt_algorithm,jet_radius);
255 }
256 break;
257 }
258
259 // Loop over all events
260 Long64_t entry, allEntries = treeReader->GetEntries();
261 cout << "** Chain contains " << allEntries << " events" << endl;
262 for(entry = 0; entry < allEntries; ++entry)
263 {
264 TLorentzVector PTmis(0,0,0,0);
265 treeReader->ReadEntry(entry);
266 treeWriter->Clear();
267 if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl;
268
269 electron.clear();
270 muon.clear();
271 elecPID.clear();
272 muonPID.clear();
273 NFCentralQ.Clear();
274
275 itGen.Reset();
276 TrackCentral.clear();
277 towers.clear();
278 input_particles.clear();
279 inclusive_jets.clear();
280 sorted_jets.clear();
281
282 // Loop over all particles in event
283 while( (particle = (TRootGenParticle*) itGen.Next()) )
284 {
285 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
286
287 int pid = abs(particle->PID);
288 float eta = fabs(particle->Eta);
289
290 //subarray of the quarks (i.e. not final particles, i.e status not equal to 1)
291 // in the tracker (i.e. for those we know the tracks)
292 // with enough p_T
293 //// This subarray is needed for the B-jet algorithm
294 // optimization for speed : put first PID condition, then ETA condition, then either pt or status
295 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 ?
296 eta < DET->MAX_TRACKER &&
297 particle->Status != 1 &&
298 particle->PT > DET->PT_QUARKS_MIN ) {
299 NFCentralQ.Add(particle);
300 }
301
302
303 // keeps only final particles, visible by the central detector, including the fiducial volume
304 // the ordering of conditions have been optimised for speed : put first the STATUS condition
305 if( (particle->Status == 1) &&
306 (
307 (pid == pMU && eta < DET->MAX_MU) ||
308 (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD)
309 )
310 ) {
311 switch(pid) {
312
313 case pE: // all electrons with eta < DET->MAX_CALO_FWD
314 DET->SmearElectron(genMomentum);
315 electron.push_back(genMomentum);
316 elecPID.push_back(particle->PID);
317 break; // case pE
318 case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD
319 DET->SmearElectron(genMomentum);
320 if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) {
321 elementPhoton = (TRootPhoton*) branchPhoton->NewEntry();
322 elementPhoton->Set(genMomentum);
323 }
324 break; // case pGAMMA
325 case pMU: // all muons with eta < DET->MAX_MU
326 DET->SmearMu(genMomentum);
327 muonPID.push_back(particle->PID);
328 muon.push_back(genMomentum);
329 break; // case pMU
330 case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD
331 case pK0S: // all K0s with eta < DET->MAX_CALO_FWD
332 DET->SmearHadron(genMomentum, 0.7);
333 break; // case hadron
334 default: // all other final particles with eta < DET->MAX_CALO_FWD
335 DET->SmearHadron(genMomentum, 1.0);
336 break;
337 } // switch (pid)
338
339 // all final particles but muons and neutrinos
340 // for calorimetric towers and mission PT
341
342 if(genMomentum.E() !=0) {
343 if(pid !=pMU) {
344 PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
345 towers.push_back(CaloTower);
346 // create a fastjet::PseudoJet with these components and put it onto
347 // back of the input_particles vector
348 input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E()));
349 //genMomentumCalo.SetPtEtaPhiE(CaloTower.Et(),CaloTower.iEta(),CaloTower.iPhi(),CaloTower.fourVector.E);
350 genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E);
351 elementCalo = (TRootCalo*) branchCalo->NewEntry();
352 elementCalo->Set(genMomentumCalo);
353 }
354 }
355
356
357 // all final charged particles
358 if(
359 ((rand()%100) < DET->TRACKING_EFF) &&
360 (genMomentum.E()!=0) &&
361 (fabs(particle->Eta) < DET->MAX_TRACKER) &&
362 (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account
363 (pid != pGAMMA) &&
364 (pid != pPI0) &&
365 (pid != pK0L) &&
366 (pid != pN) &&
367 (pid != pSIGMA0) &&
368 (pid != pDELTA0) &&
369 (pid != pK0S) // not charged particles : invisible by tracker
370 )
371 {
372 elementTracks = (TRootTracks*) branchTracks->NewEntry();
373 elementTracks->Set(genMomentum);
374 TrackCentral.push_back(genMomentum);
375 }
376
377} // switch
378
379 // Forward particles in CASTOR ?
380 /* if (particle->Status == 1 && (fabs(particle->Eta) > DET->MIN_CALO_VFWD)
381 && (fabs(particle->Eta) < DET->MAX_CALO_VFWD)) {
382
383
384 } // CASTOR
385 */
386 // Zero degree calorimeter, for forward neutrons and photons
387 if (particle->Status ==1 && (pid == pN || pid == pGAMMA ) && eta > DET->MIN_ZDC ) {
388 // !!!!!!!!! vérifier que particle->Z est bien en micromÚtres!!!
389 // !!!!!!!!! vérifier que particle->T est bien en secondes!!!
390 // !!!!!!!!! pas de smearing ! on garde trop d'info !
391 elementZdc = (TRootZdcHits*) branchZDC->NewEntry();
392 elementZdc->Set(genMomentum);
393
394 // time of flight t is t = T + d/[ cos(theta) v ]
395 //double tx = acos(particle->Px/particle->Pz);
396 //double ty = acos(particle->Py/particle->Pz);
397 //double theta = (1E-6)*sqrt( pow(tx,2) + pow(ty,2) );
398 //double flight_distance = (DET->ZDC_S - particle->Z*(1E-6))/cos(theta) ; // assumes that Z is in micrometers
399 double flight_distance = (DET->ZDC_S - particle->Z*(1E-6));
400 // assumes also that the emission angle is so small that 1/(cos theta) = 1
401 elementZdc->T = 0*particle->T + flight_distance/speed_of_light; // assumes highly relativistic particles
402 elementZdc->side = sign(particle->Eta);
403
404 }
405
406 // if forward proton
407 if( (pid == pP) && (particle->Status == 1) && (fabs(particle->Eta) > DET->MAX_CALO_FWD) )
408 {
409 // !!!!!!!! put here particle->CHARGE and particle->MASS
410 H_BeamParticle p1; /// put here particle->CHARGE and particle->MASS
411 p1.smearAng();
412 p1.smearPos();
413 p1.setPosition(p1.getX()-500.,p1.getY(),p1.getTX()-1*kickers_on*CRANG,p1.getTY(),0);
414 p1.set4Momentum(particle->Px,particle->Py,particle->Pz,particle->E);
415
416 H_BeamLine *beamline;
417 if(particle->Eta >0) beamline = beamline1;
418 else beamline = beamline2;
419
420 p1.computePath(beamline,1);
421
422 if(p1.stopped(beamline)) {
423 if (p1.getStoppingElement()->getName()=="rp220_1" || p1.getStoppingElement()->getName()=="rp220_2") {
424 p1.propagate(DET->RP220_S);
425 elementRP220 = (TRootRomanPotHits*) branchRP220->NewEntry();
426 elementRP220->X = (1E-6)*p1.getX(); // [m]
427 elementRP220->Y = (1E-6)*p1.getY(); // [m]
428 elementRP220->Tx = (1E-6)*p1.getTX(); // [rad]
429 elementRP220->Ty = (1E-6)*p1.getTY(); // [rad]
430 elementRP220->S = p1.getS(); // [m]
431 elementRP220->T = -1; // not yet implemented
432 elementRP220->E = p1.getE(); // not yet implemented
433 elementRP220->q2 = -1; // not yet implemented
434 elementRP220->side = sign(particle->Eta);
435
436 } else if (p1.getStoppingElement()->getName()=="rp420_1" || p1.getStoppingElement()->getName()=="rp420_2") {
437 p1.propagate(DET->FP420_S);
438 elementFP420 = (TRootRomanPotHits*) branchFP420->NewEntry();
439 elementFP420->X = (1E-6)*p1.getX(); // [m]
440 elementFP420->Y = (1E-6)*p1.getY(); // [m]
441 elementFP420->Tx = (1E-6)*p1.getTX(); // [rad]
442 elementFP420->Ty = (1E-6)*p1.getTY(); // [rad]
443 elementFP420->S = p1.getS(); // [m]
444 elementFP420->T = -1; // not yet implemented
445 elementFP420->E = p1.getE(); // not yet implemented
446 elementFP420->q2 = -1; // not yet implemented
447 elementFP420->side = sign(particle->Eta);
448 }
449 }
450
451 // if(p1.stopped(beamline) && (p1.getStoppingElement()->getS() > 100))
452 // cout << "Eloss =" << 7000.-p1.getE() << " ; " << p1.getStoppingElement()->getName() << endl;
453 } // if forward proton
454
455 } // while
456
457 for(unsigned int i=0; i < electron.size(); i++) {
458 if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt)
459 {
460 elementElec = (TRootElectron*) branchElectron->NewEntry();
461 elementElec->Set(electron[i]);
462 elementElec->Charge = sign(elecPID[i]);
463 elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0);
464 }
465 }
466 for(unsigned int i=0; i < muon.size(); i++) {
467 if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt)
468 {
469 elementMu = (TRootMuon*) branchMuon->NewEntry();
470 elementMu->Charge = sign(muonPID[i]);
471 elementMu->Set(muon[i]);
472 elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0);
473 }
474 }
475
476 // computes the Missing Transverse Momentum
477 TLorentzVector Att(0.,0.,0.,0.);
478 for(unsigned int i=0; i < towers.size(); i++)
479 {
480 Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E);
481 PTmis = PTmis + Att;
482 }
483 elementEtmis = (TRootETmis*) branchETmis->NewEntry();
484 elementEtmis->ET = (PTmis).Pt();
485 elementEtmis->Phi = (-PTmis).Phi();
486 elementEtmis->Px = (-PTmis).Px();
487 elementEtmis->Py = (-PTmis).Py();
488 //*****************************
489
490 // run the jet clustering with the above jet definition
491 if(input_particles.size()!=0)
492 {
493 fastjet::ClusterSequence clust_seq(input_particles, jet_def);
494 // extract the inclusive jets with pt > 5 GeV
495 double ptmin = 5.0;
496 inclusive_jets = clust_seq.inclusive_jets(ptmin);
497 // sort jets into increasing pt
498 sorted_jets = sorted_by_pt(inclusive_jets);
499 }
500 for (unsigned int i = 0; i < sorted_jets.size(); i++) {
501 TLorentzVector JET;
502 JET.SetPxPyPzE(sorted_jets[i].px(),sorted_jets[i].py(),sorted_jets[i].pz(),sorted_jets[i].E());
503 // Tau jet identification : 1! track and electromagnetic collimation
504 if(fabs(JET.Eta()) < (DET->MAX_TRACKER - DET->TAU_CONE_TRACKS)) {
505 double Energie_tau_central = DET->EnergySmallCone(towers,JET.Eta(),JET.Phi());
506 if(
507 ( Energie_tau_central/JET.E() > DET->TAU_EM_COLLIMATION ) &&
508 ( DET->NumTracks(TrackCentral,DET->PT_TRACK_TAU,JET.Eta(),JET.Phi()) == 1 ) &&
509 ( JET.Pt() > DET->TAUJET_pt)
510 ) {
511 elementTauJet = (TRootTauJet*) branchTauJet->NewEntry();
512 elementTauJet->Set(JET);
513 } // if tau jet
514 } // if JET.eta < tracker - tau_cone : Tau jet identification
515
516 if(JET.Pt() > DET->JET_pt)
517 {
518 elementJet = (TRootJet*) branchJet->NewEntry();
519 elementJet->Set(JET);
520 // b-jets
521 bool btag=false;
522 if((fabs(JET.Eta()) < DET->MAX_TRACKER && DET->Btaggedjet(JET, NFCentralQ)))btag=true;
523 elementJet->Btag = btag;
524 }
525 } // for itJet : loop on all jets
526
527 treeWriter->Fill();
528 // Add here the trigger
529 // Should test all the trigger table on the event, based on reconstructed objects
530 } // Loop over all events
531 treeWriter->Write();
532
533 cout << "** Exiting..." << endl;
534
535 delete treeWriter;
536 delete treeReader;
537 delete DET;
538 if(converter) delete converter;
539
540 todo("TODO");
541}
542
Note: See TracBrowser for help on using the repository browser.