Fork me on GitHub

source: git/modules/TreeWriter.cc@ b0fa9af

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b0fa9af was acd0621, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

including tracking variables

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