Fork me on GitHub

source: svn/trunk/modules/TreeWriter.cc

Last change on this file was 1372, checked in by Pavel Demin, 10 years ago

several minor corrections

  • Property svn:keywords set to Id Revision Date
File size: 18.8 KB
RevLine 
[694]1
[815]2/** \class TreeWriter
3 *
4 * Fills ROOT tree branches.
5 *
6 * $Date: 2014-04-17 08:54:17 +0000 (Thu, 17 Apr 2014) $
7 * $Revision: 1372 $
8 *
9 *
10 * \author P. Demin - UCL, Louvain-la-Neuve
11 *
12 */
13
[694]14#include "modules/TreeWriter.h"
15
[703]16#include "classes/DelphesClasses.h"
17#include "classes/DelphesFactory.h"
[766]18#include "classes/DelphesFormula.h"
[694]19
20#include "ExRootAnalysis/ExRootResult.h"
[703]21#include "ExRootAnalysis/ExRootFilter.h"
22#include "ExRootAnalysis/ExRootClassifier.h"
[694]23#include "ExRootAnalysis/ExRootTreeBranch.h"
24
25#include "TROOT.h"
[703]26#include "TMath.h"
[694]27#include "TString.h"
[703]28#include "TFormula.h"
29#include "TRandom3.h"
30#include "TObjArray.h"
31#include "TDatabasePDG.h"
[694]32#include "TLorentzVector.h"
33
[988]34#include <algorithm>
[703]35#include <stdexcept>
[694]36#include <iostream>
[703]37#include <sstream>
[694]38
39using namespace std;
40
41//------------------------------------------------------------------------------
42
43TreeWriter::TreeWriter()
44{
45}
46
47//------------------------------------------------------------------------------
48
49TreeWriter::~TreeWriter()
50{
51}
52
53//------------------------------------------------------------------------------
54
55void TreeWriter::Init()
56{
57 fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
[1325]58 fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
[694]59 fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
60 fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
61 fClassMap[Photon::Class()] = &TreeWriter::ProcessPhotons;
62 fClassMap[Electron::Class()] = &TreeWriter::ProcessElectrons;
63 fClassMap[Muon::Class()] = &TreeWriter::ProcessMuons;
64 fClassMap[Jet::Class()] = &TreeWriter::ProcessJets;
65 fClassMap[MissingET::Class()] = &TreeWriter::ProcessMissingET;
[893]66 fClassMap[ScalarHT::Class()] = &TreeWriter::ProcessScalarHT;
[1114]67 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
[1123]68 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
[1361]69 fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
[694]70
71 TBranchMap::iterator itBranchMap;
72 map< TClass *, TProcessMethod >::iterator itClassMap;
73
74 // read branch configuration and
75 // import array with output from filter/classifier/jetfinder modules
76
77 ExRootConfParam param = GetParam("Branch");
78 Long_t i, size;
79 TString branchName, branchClassName, branchInputArray;
80 TClass *branchClass;
[895]81 TObjArray *array;
[694]82 ExRootTreeBranch *branch;
83
84 size = param.GetSize();
[703]85 for(i = 0; i < size/3; ++i)
[694]86 {
[703]87 branchInputArray = param[i*3].GetString();
88 branchName = param[i*3 + 1].GetString();
89 branchClassName = param[i*3 + 2].GetString();
[694]90
91 branchClass = gROOT->GetClass(branchClassName);
92
93 if(!branchClass)
94 {
95 cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
96 continue;
97 }
[1361]98
[694]99 itClassMap = fClassMap.find(branchClass);
100 if(itClassMap == fClassMap.end())
101 {
102 cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
103 continue;
104 }
105
106 array = ImportArray(branchInputArray);
107 branch = NewBranch(branchName, branchClass);
108
[895]109 fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
[694]110 }
111
112}
113
114//------------------------------------------------------------------------------
115
116void TreeWriter::Finish()
117{
118}
119
120//------------------------------------------------------------------------------
121
[930]122void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
123{
124 TIter it1(candidate->GetCandidates());
125 it1.Reset();
126 array->Clear();
127 while((candidate = static_cast<Candidate*>(it1.Next())))
128 {
129 TIter it2(candidate->GetCandidates());
130
131 // particle
132 if(candidate->GetCandidates()->GetEntriesFast() == 0)
133 {
134 array->Add(candidate);
135 continue;
136 }
137
138 // track
139 candidate = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
140 if(candidate->GetCandidates()->GetEntriesFast() == 0)
141 {
142 array->Add(candidate);
143 continue;
144 }
145
146 // tower
147 it2.Reset();
148 while((candidate = static_cast<Candidate*>(it2.Next())))
149 {
150 array->Add(candidate->GetCandidates()->At(0));
151 }
152 }
153}
154
155//------------------------------------------------------------------------------
156
[895]157void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
[694]158{
[895]159 TIter iterator(array);
[694]160 Candidate *candidate = 0;
161 GenParticle *entry = 0;
[771]162 Double_t pt, signPz, cosTheta, eta, rapidity;
[1361]163
[1345]164 const Double_t c_light = 2.99792458E8;
[1361]165
[694]166 // loop over all particles
[895]167 iterator.Reset();
168 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]169 {
170 const TLorentzVector &momentum = candidate->Momentum;
171 const TLorentzVector &position = candidate->Position;
172
173 entry = static_cast<GenParticle*>(branch->NewEntry());
174
[921]175 entry->SetBit(kIsReferenced);
176 entry->SetUniqueID(candidate->GetUniqueID());
[923]177
[694]178 pt = momentum.Pt();
[771]179 cosTheta = TMath::Abs(momentum.CosTheta());
[694]180 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]181 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
182 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]183
184 entry->PID = candidate->PID;
185
186 entry->Status = candidate->Status;
[1014]187 entry->IsPU = candidate->IsPU;
[694]188
189 entry->M1 = candidate->M1;
190 entry->M2 = candidate->M2;
191
[845]192 entry->D1 = candidate->D1;
193 entry->D2 = candidate->D2;
194
[694]195 entry->Charge = candidate->Charge;
196 entry->Mass = candidate->Mass;
197
198 entry->E = momentum.E();
199 entry->Px = momentum.Px();
200 entry->Py = momentum.Py();
201 entry->Pz = momentum.Pz();
202
203 entry->Eta = eta;
204 entry->Phi = momentum.Phi();
205 entry->PT = pt;
206
207 entry->Rapidity = rapidity;
208
209 entry->X = position.X();
210 entry->Y = position.Y();
211 entry->Z = position.Z();
[1345]212 entry->T = position.T()*1.0E-3/c_light;
[694]213 }
214}
215
216//------------------------------------------------------------------------------
217
[1323]218void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
219{
220 TIter iterator(array);
221 Candidate *candidate = 0;
222 Vertex *entry = 0;
[1361]223
[1345]224 const Double_t c_light = 2.99792458E8;
[1323]225
226 // loop over all vertices
227 iterator.Reset();
228 while((candidate = static_cast<Candidate*>(iterator.Next())))
229 {
230 const TLorentzVector &position = candidate->Position;
231
232 entry = static_cast<Vertex*>(branch->NewEntry());
233
234 entry->X = position.X();
235 entry->Y = position.Y();
236 entry->Z = position.Z();
[1345]237 entry->T = position.T()*1.0E-3/c_light;
[1323]238 }
239}
240
241//------------------------------------------------------------------------------
242
[895]243void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
[694]244{
[895]245 TIter iterator(array);
[694]246 Candidate *candidate = 0;
[920]247 Candidate *particle = 0;
[694]248 Track *entry = 0;
[771]249 Double_t pt, signz, cosTheta, eta, rapidity;
[1345]250 const Double_t c_light = 2.99792458E8;
[1361]251
[1323]252 // loop over all tracks
[895]253 iterator.Reset();
254 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]255 {
[920]256 const TLorentzVector &position = candidate->Position;
[694]257
[920]258 cosTheta = TMath::Abs(position.CosTheta());
259 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
260 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
261 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
[694]262
263 entry = static_cast<Track*>(branch->NewEntry());
264
[936]265 entry->SetBit(kIsReferenced);
266 entry->SetUniqueID(candidate->GetUniqueID());
267
[703]268 entry->PID = candidate->PID;
269
270 entry->Charge = candidate->Charge;
271
[694]272 entry->EtaOuter = eta;
[920]273 entry->PhiOuter = position.Phi();
[696]274
[920]275 entry->XOuter = position.X();
276 entry->YOuter = position.Y();
277 entry->ZOuter = position.Z();
[1345]278 entry->TOuter = position.T()*1.0E-3/c_light;
[1367]279
280 entry->Dxy = candidate->Dxy;
281 entry->SDxy = candidate->SDxy ;
[1372]282 entry->Xd = candidate->Xd;
283 entry->Yd = candidate->Yd;
284 entry->Zd = candidate->Zd;
[1367]285
[920]286 const TLorentzVector &momentum = candidate->Momentum;
[694]287
288 pt = momentum.Pt();
[771]289 cosTheta = TMath::Abs(momentum.CosTheta());
[694]290 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]291 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
292 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
[694]293
294 entry->Eta = eta;
295 entry->Phi = momentum.Phi();
296 entry->PT = pt;
[988]297
[920]298 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
299 const TLorentzVector &initialPosition = particle->Position;
300
301 entry->X = initialPosition.X();
302 entry->Y = initialPosition.Y();
303 entry->Z = initialPosition.Z();
[1345]304 entry->T = initialPosition.T()*1.0E-3/c_light;
[988]305
[920]306 entry->Particle = particle;
[694]307 }
308}
309
310//------------------------------------------------------------------------------
311
[895]312void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
[694]313{
[895]314 TIter iterator(array);
[694]315 Candidate *candidate = 0;
316 Tower *entry = 0;
[771]317 Double_t pt, signPz, cosTheta, eta, rapidity;
[1345]318 const Double_t c_light = 2.99792458E8;
[1361]319
[1323]320 // loop over all towers
[895]321 iterator.Reset();
322 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]323 {
324 const TLorentzVector &momentum = candidate->Momentum;
[1345]325 const TLorentzVector &position = candidate->Position;
[1361]326
[694]327 pt = momentum.Pt();
[771]328 cosTheta = TMath::Abs(momentum.CosTheta());
[694]329 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]330 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
331 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]332
333 entry = static_cast<Tower*>(branch->NewEntry());
334
[936]335 entry->SetBit(kIsReferenced);
336 entry->SetUniqueID(candidate->GetUniqueID());
337
[694]338 entry->Eta = eta;
339 entry->Phi = momentum.Phi();
340 entry->ET = pt;
[696]341 entry->E = momentum.E();
342 entry->Eem = candidate->Eem;
343 entry->Ehad = candidate->Ehad;
[898]344 entry->Edges[0] = candidate->Edges[0];
345 entry->Edges[1] = candidate->Edges[1];
346 entry->Edges[2] = candidate->Edges[2];
347 entry->Edges[3] = candidate->Edges[3];
[1361]348
[1345]349 entry->T = position.T()*1.0E-3/c_light;
[1361]350
[930]351 FillParticles(candidate, &entry->Particles);
[694]352 }
353}
354
355//------------------------------------------------------------------------------
356
[895]357void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
[694]358{
[895]359 TIter iterator(array);
[694]360 Candidate *candidate = 0;
361 Photon *entry = 0;
[771]362 Double_t pt, signPz, cosTheta, eta, rapidity;
[1345]363 const Double_t c_light = 2.99792458E8;
[1361]364
[1066]365 array->Sort();
366
367 // loop over all photons
[895]368 iterator.Reset();
369 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]370 {
[930]371 TIter it1(candidate->GetCandidates());
[694]372 const TLorentzVector &momentum = candidate->Momentum;
[1345]373 const TLorentzVector &position = candidate->Position;
[694]374
[1361]375
[694]376 pt = momentum.Pt();
[771]377 cosTheta = TMath::Abs(momentum.CosTheta());
[694]378 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]379 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
380 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]381
382 entry = static_cast<Photon*>(branch->NewEntry());
383
384 entry->Eta = eta;
385 entry->Phi = momentum.Phi();
386 entry->PT = pt;
[923]387 entry->E = momentum.E();
[1361]388
[1345]389 entry->T = position.T()*1.0E-3/c_light;
[1361]390
[988]391 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
392
[930]393 FillParticles(candidate, &entry->Particles);
[694]394 }
395}
396
397//------------------------------------------------------------------------------
398
[895]399void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
[694]400{
[895]401 TIter iterator(array);
[694]402 Candidate *candidate = 0;
403 Electron *entry = 0;
[771]404 Double_t pt, signPz, cosTheta, eta, rapidity;
[1345]405 const Double_t c_light = 2.99792458E8;
[1361]406
[1066]407 array->Sort();
408
409 // loop over all electrons
[895]410 iterator.Reset();
411 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]412 {
413 const TLorentzVector &momentum = candidate->Momentum;
[1345]414 const TLorentzVector &position = candidate->Position;
[1361]415
[694]416 pt = momentum.Pt();
[771]417 cosTheta = TMath::Abs(momentum.CosTheta());
[694]418 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]419 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
420 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]421
422 entry = static_cast<Electron*>(branch->NewEntry());
423
424 entry->Eta = eta;
425 entry->Phi = momentum.Phi();
426 entry->PT = pt;
[1361]427
[1345]428 entry->T = position.T()*1.0E-3/c_light;
[1361]429
[694]430 entry->Charge = candidate->Charge;
[920]431
[988]432 entry->EhadOverEem = 0.0;
433
[920]434 entry->Particle = candidate->GetCandidates()->At(0);
[694]435 }
436}
437
438//------------------------------------------------------------------------------
439
[895]440void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
[694]441{
[895]442 TIter iterator(array);
[694]443 Candidate *candidate = 0;
444 Muon *entry = 0;
[771]445 Double_t pt, signPz, cosTheta, eta, rapidity;
[1361]446
[1345]447 const Double_t c_light = 2.99792458E8;
[1361]448
[1066]449 array->Sort();
450
451 // loop over all muons
[895]452 iterator.Reset();
453 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]454 {
455 const TLorentzVector &momentum = candidate->Momentum;
[1345]456 const TLorentzVector &position = candidate->Position;
[694]457
[1361]458
[694]459 pt = momentum.Pt();
[771]460 cosTheta = TMath::Abs(momentum.CosTheta());
[694]461 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]462 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
463 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]464
465 entry = static_cast<Muon*>(branch->NewEntry());
466
[939]467 entry->SetBit(kIsReferenced);
468 entry->SetUniqueID(candidate->GetUniqueID());
469
[694]470 entry->Eta = eta;
471 entry->Phi = momentum.Phi();
472 entry->PT = pt;
473
[1345]474 entry->T = position.T()*1.0E-3/c_light;
[1361]475
[694]476 entry->Charge = candidate->Charge;
[920]477
478 entry->Particle = candidate->GetCandidates()->At(0);
[694]479 }
480}
481
482//------------------------------------------------------------------------------
483
[895]484void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
[694]485{
[895]486 TIter iterator(array);
[936]487 Candidate *candidate = 0, *constituent = 0;
[694]488 Jet *entry = 0;
[771]489 Double_t pt, signPz, cosTheta, eta, rapidity;
[988]490 Double_t ecalEnergy, hcalEnergy;
[1345]491 const Double_t c_light = 2.99792458E8;
[1361]492
[1066]493 array->Sort();
494
[694]495 // loop over all jets
[895]496 iterator.Reset();
497 while((candidate = static_cast<Candidate*>(iterator.Next())))
[694]498 {
[936]499 TIter itConstituents(candidate->GetCandidates());
[1361]500
[694]501 const TLorentzVector &momentum = candidate->Momentum;
[1345]502 const TLorentzVector &position = candidate->Position;
[1361]503
[694]504 pt = momentum.Pt();
[771]505 cosTheta = TMath::Abs(momentum.CosTheta());
[694]506 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[771]507 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
508 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
[694]509
510 entry = static_cast<Jet*>(branch->NewEntry());
511
512 entry->Eta = eta;
513 entry->Phi = momentum.Phi();
514 entry->PT = pt;
515
[1345]516 entry->T = position.T()*1.0E-3/c_light;
[1361]517
[694]518 entry->Mass = momentum.M();
[696]519
[900]520 entry->DeltaEta = candidate->DeltaEta;
521 entry->DeltaPhi = candidate->DeltaPhi;
[696]522
523 entry->BTag = candidate->BTag;
[926]524 entry->TauTag = candidate->TauTag;
[900]525
[926]526 entry->Charge = candidate->Charge;
527
[936]528 itConstituents.Reset();
529 entry->Constituents.Clear();
[988]530 ecalEnergy = 0.0;
531 hcalEnergy = 0.0;
[936]532 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
533 {
534 entry->Constituents.Add(constituent);
[988]535 ecalEnergy += constituent->Eem;
536 hcalEnergy += constituent->Ehad;
[936]537 }
[988]538
539 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
[1361]540
[1348]541 //--- Pile-Up Jet ID variables ----
[988]542
[1348]543 entry->NCharged = candidate->NCharged;
544 entry->NNeutrals = candidate->NNeutrals;
545 entry->Beta = candidate->Beta;
546 entry->BetaStar = candidate->BetaStar;
547 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
548 entry->PTD = candidate->PTD;
549 entry->FracPt[0] = candidate->FracPt[0];
550 entry->FracPt[1] = candidate->FracPt[1];
551 entry->FracPt[2] = candidate->FracPt[2];
552 entry->FracPt[3] = candidate->FracPt[3];
553 entry->FracPt[4] = candidate->FracPt[4];
[1368]554
555 //--- N-subjettiness variables ----
556
[1372]557 entry->Tau1 = candidate->Tau[0];
558 entry->Tau2 = candidate->Tau[1];
559 entry->Tau3 = candidate->Tau[2];
560 entry->Tau4 = candidate->Tau[3];
561 entry->Tau5 = candidate->Tau[4];
[1368]562
[930]563 FillParticles(candidate, &entry->Particles);
[694]564 }
565}
566
567//------------------------------------------------------------------------------
568
[895]569void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
[694]570{
571 Candidate *candidate = 0;
572 MissingET *entry = 0;
573
[893]574 // get the first entry
[895]575 if((candidate = static_cast<Candidate*>(array->At(0))))
[694]576 {
577 const TLorentzVector &momentum = candidate->Momentum;
578
579 entry = static_cast<MissingET*>(branch->NewEntry());
580
[1321]581 entry->Eta = (-momentum).Eta();
[1092]582 entry->Phi = (-momentum).Phi();
[893]583 entry->MET = momentum.Pt();
[694]584 }
585}
586
587//------------------------------------------------------------------------------
588
[895]589void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
[893]590{
591 Candidate *candidate = 0;
592 ScalarHT *entry = 0;
593
594 // get the first entry
[895]595 if((candidate = static_cast<Candidate*>(array->At(0))))
[893]596 {
597 const TLorentzVector &momentum = candidate->Momentum;
598
599 entry = static_cast<ScalarHT*>(branch->NewEntry());
600
601 entry->HT = momentum.Pt();
602 }
603}
604
605//------------------------------------------------------------------------------
606
[1114]607void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
608{
[1321]609 TIter iterator(array);
[1114]610 Candidate *candidate = 0;
611 Rho *entry = 0;
612
[1321]613 // loop over all rho
614 iterator.Reset();
615 while((candidate = static_cast<Candidate*>(iterator.Next())))
[1114]616 {
617 const TLorentzVector &momentum = candidate->Momentum;
618
619 entry = static_cast<Rho*>(branch->NewEntry());
620
[1115]621 entry->Rho = momentum.E();
[1321]622 entry->Edges[0] = candidate->Edges[0];
623 entry->Edges[1] = candidate->Edges[1];
[1114]624 }
625}
626
627//------------------------------------------------------------------------------
628
[1123]629void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
630{
631 Candidate *candidate = 0;
632 Weight *entry = 0;
633
634 // get the first entry
635 if((candidate = static_cast<Candidate*>(array->At(0))))
636 {
637 const TLorentzVector &momentum = candidate->Momentum;
638
639 entry = static_cast<Weight*>(branch->NewEntry());
640
641 entry->Weight = momentum.E();
642 }
643}
644
645//------------------------------------------------------------------------------
646
[1361]647void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
648{
649 TIter iterator(array);
650 Candidate *candidate = 0;
651 HectorHit *entry = 0;
652
653 // loop over all roman pot hits
654 iterator.Reset();
655 while((candidate = static_cast<Candidate*>(iterator.Next())))
656 {
657 const TLorentzVector &position = candidate->Position;
658 const TLorentzVector &momentum = candidate->Momentum;
659
660 entry = static_cast<HectorHit*>(branch->NewEntry());
661
662 entry->E = momentum.E();
663
664 entry->Tx = momentum.Px();
665 entry->Ty = momentum.Py();
666
667 entry->T = position.T();
668
669 entry->X = position.X();
670 entry->Y = position.Y();
671 entry->S = position.Z();
[1366]672
673 entry->Particle = candidate->GetCandidates()->At(0);
[1361]674 }
675}
676
677//------------------------------------------------------------------------------
678
[694]679void TreeWriter::Process()
680{
681 TBranchMap::iterator itBranchMap;
682 ExRootTreeBranch *branch;
683 TProcessMethod method;
[895]684 TObjArray *array;
[694]685
686 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
687 {
688 branch = itBranchMap->first;
689 method = itBranchMap->second.first;
[895]690 array = itBranchMap->second.second;
[694]691
[895]692 (this->*method)(branch, array);
[694]693 }
694}
695
696//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.