Fork me on GitHub

source: svn/trunk/Delphes.cpp@ 38

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

remove MY bug jets

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