Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 53

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

2 new jet algorithms

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