Fork me on GitHub

source: git/modules/TreeWriter.cc@ 5496767

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

added array of constituents to vertices

  • Property mode set to 100644
File size: 22.7 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, *constituent = 0;
243 Vertex *entry = 0;
244
245 const Double_t c_light = 2.99792458E8;
246
247 Double_t x, y, z, t, xError, yError, zError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2;
248 UInt_t index, ndf;
249
250 array->Sort();
251
252 // loop over all vertices
253 iterator.Reset();
254 while((candidate = static_cast<Candidate*>(iterator.Next())))
255 {
256
257 index = candidate->ClusterIndex;
258 ndf = candidate->ClusterNDF;
259 sigma = candidate->ClusterSigma;
260 sumPT2 = candidate->SumPT2;
261 btvSumPT2 = candidate->BTVSumPT2;
262 genDeltaZ = candidate->GenDeltaZ;
263 genSumPT2 = candidate->GenSumPT2;
264
265 x = candidate->Position.X();
266 y = candidate->Position.Y();
267 z = candidate->Position.Z();
268 t = candidate->Position.T()*1.0E-3/c_light;
269
270 xError = candidate->PositionError.X ();
271 yError = candidate->PositionError.Y ();
272 zError = candidate->PositionError.Z ();
273
274 entry = static_cast<Vertex*>(branch->NewEntry());
275
276 entry->Index = index;
277 entry->NDF = ndf;
278 entry->Sigma = sigma;
279 entry->SumPT2 = sumPT2;
280 entry->BTVSumPT2 = btvSumPT2;
281 entry->GenDeltaZ = genDeltaZ;
282 entry->GenSumPT2 = genSumPT2;
283
284 entry->X = x;
285 entry->Y = y;
286 entry->Z = z;
287 entry->T = t;
288
289 entry->ErrorX = xError;
290 entry->ErrorY = yError;
291 entry->ErrorZ = zError;
292
293 TIter itConstituents(candidate->GetCandidates());
294
295 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
296 {
297 entry->Constituents.Add(constituent);
298 }
299
300 }
301}
302
303
304//------------------------------------------------------------------------------
305
306void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
307{
308 TIter iterator(array);
309 Candidate *candidate = 0;
310 Candidate *particle = 0;
311 Track *entry = 0;
312 Double_t pt, signz, cosTheta, eta, rapidity;
313 const Double_t c_light = 2.99792458E8;
314
315 // loop over all tracks
316 iterator.Reset();
317 while((candidate = static_cast<Candidate*>(iterator.Next())))
318 {
319 const TLorentzVector &position = candidate->Position;
320
321 cosTheta = TMath::Abs(position.CosTheta());
322 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
323 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
324 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
325
326 entry = static_cast<Track*>(branch->NewEntry());
327
328 entry->SetBit(kIsReferenced);
329 entry->SetUniqueID(candidate->GetUniqueID());
330
331 entry->PID = candidate->PID;
332
333 entry->Charge = candidate->Charge;
334
335 entry->EtaOuter = eta;
336 entry->PhiOuter = position.Phi();
337
338 entry->XOuter = position.X();
339 entry->YOuter = position.Y();
340 entry->ZOuter = position.Z();
341 entry->TOuter = position.T()*1.0E-3/c_light;
342
343 entry->L = candidate->L;
344
345 entry->D0 = candidate->D0;
346 entry->ErrorD0 = candidate->ErrorD0;
347 entry->DZ = candidate->DZ;
348 entry->ErrorDZ = candidate->ErrorDZ;
349 entry->P = candidate->P;
350 entry->ErrorP = candidate->ErrorP;
351 entry->PT = candidate->PT;
352 entry->ErrorPT = candidate->ErrorPT;
353 entry->CtgTheta = candidate->CtgTheta;
354 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
355 entry->Phi = candidate->Phi;
356 entry->ErrorPhi = candidate->ErrorPhi;
357
358 entry->Xd = candidate->Xd;
359 entry->Yd = candidate->Yd;
360 entry->Zd = candidate->Zd;
361
362 const TLorentzVector &momentum = candidate->Momentum;
363
364 pt = momentum.Pt();
365 cosTheta = TMath::Abs(momentum.CosTheta());
366 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
367 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
368 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
369
370 entry->Eta = eta;
371
372 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
373 const TLorentzVector &initialPosition = particle->Position;
374
375 entry->X = initialPosition.X();
376 entry->Y = initialPosition.Y();
377 entry->Z = initialPosition.Z();
378 entry->T = initialPosition.T()*1.0E-3/c_light;
379
380 entry->Particle = particle;
381
382 entry->VertexIndex = candidate->ClusterIndex;
383
384 }
385}
386
387//------------------------------------------------------------------------------
388
389void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
390{
391 TIter iterator(array);
392 Candidate *candidate = 0;
393 Tower *entry = 0;
394 Double_t pt, signPz, cosTheta, eta, rapidity;
395 const Double_t c_light = 2.99792458E8;
396
397 // loop over all towers
398 iterator.Reset();
399 while((candidate = static_cast<Candidate*>(iterator.Next())))
400 {
401 const TLorentzVector &momentum = candidate->Momentum;
402 const TLorentzVector &position = candidate->Position;
403
404 pt = momentum.Pt();
405 cosTheta = TMath::Abs(momentum.CosTheta());
406 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
407 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
408 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
409
410 entry = static_cast<Tower*>(branch->NewEntry());
411
412 entry->SetBit(kIsReferenced);
413 entry->SetUniqueID(candidate->GetUniqueID());
414
415 entry->Eta = eta;
416 entry->Phi = momentum.Phi();
417 entry->ET = pt;
418 entry->E = momentum.E();
419 entry->Eem = candidate->Eem;
420 entry->Ehad = candidate->Ehad;
421 entry->Edges[0] = candidate->Edges[0];
422 entry->Edges[1] = candidate->Edges[1];
423 entry->Edges[2] = candidate->Edges[2];
424 entry->Edges[3] = candidate->Edges[3];
425
426 entry->T = position.T()*1.0E-3/c_light;
427 entry->NTimeHits = candidate->NTimeHits;
428
429 FillParticles(candidate, &entry->Particles);
430 }
431}
432
433//------------------------------------------------------------------------------
434
435void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
436{
437 TIter iterator(array);
438 Candidate *candidate = 0;
439 Photon *entry = 0;
440 Double_t pt, signPz, cosTheta, eta, rapidity;
441 const Double_t c_light = 2.99792458E8;
442
443 array->Sort();
444
445 // loop over all photons
446 iterator.Reset();
447 while((candidate = static_cast<Candidate*>(iterator.Next())))
448 {
449 TIter it1(candidate->GetCandidates());
450 const TLorentzVector &momentum = candidate->Momentum;
451 const TLorentzVector &position = candidate->Position;
452
453
454 pt = momentum.Pt();
455 cosTheta = TMath::Abs(momentum.CosTheta());
456 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
457 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
458 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
459
460 entry = static_cast<Photon*>(branch->NewEntry());
461
462 entry->Eta = eta;
463 entry->Phi = momentum.Phi();
464 entry->PT = pt;
465 entry->E = momentum.E();
466
467 entry->T = position.T()*1.0E-3/c_light;
468
469 // Isolation variables
470
471 entry->IsolationVar = candidate->IsolationVar;
472 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
473 entry->SumPtCharged = candidate->SumPtCharged ;
474 entry->SumPtNeutral = candidate->SumPtNeutral ;
475 entry->SumPtChargedPU = candidate->SumPtChargedPU ;
476 entry->SumPt = candidate->SumPt ;
477
478 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
479
480 FillParticles(candidate, &entry->Particles);
481 }
482}
483
484//------------------------------------------------------------------------------
485
486void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
487{
488 TIter iterator(array);
489 Candidate *candidate = 0;
490 Electron *entry = 0;
491 Double_t pt, signPz, cosTheta, eta, rapidity;
492 const Double_t c_light = 2.99792458E8;
493
494 array->Sort();
495
496 // loop over all electrons
497 iterator.Reset();
498 while((candidate = static_cast<Candidate*>(iterator.Next())))
499 {
500 const TLorentzVector &momentum = candidate->Momentum;
501 const TLorentzVector &position = candidate->Position;
502
503 pt = momentum.Pt();
504 cosTheta = TMath::Abs(momentum.CosTheta());
505 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
506 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
507 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
508
509 entry = static_cast<Electron*>(branch->NewEntry());
510
511 entry->Eta = eta;
512 entry->Phi = momentum.Phi();
513 entry->PT = pt;
514
515 entry->T = position.T()*1.0E-3/c_light;
516
517 // Isolation variables
518
519 entry->IsolationVar = candidate->IsolationVar;
520 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
521 entry->SumPtCharged = candidate->SumPtCharged ;
522 entry->SumPtNeutral = candidate->SumPtNeutral ;
523 entry->SumPtChargedPU = candidate->SumPtChargedPU ;
524 entry->SumPt = candidate->SumPt ;
525
526
527 entry->Charge = candidate->Charge;
528
529 entry->EhadOverEem = 0.0;
530
531 entry->Particle = candidate->GetCandidates()->At(0);
532 }
533}
534
535//------------------------------------------------------------------------------
536
537void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
538{
539 TIter iterator(array);
540 Candidate *candidate = 0;
541 Muon *entry = 0;
542 Double_t pt, signPz, cosTheta, eta, rapidity;
543
544 const Double_t c_light = 2.99792458E8;
545
546 array->Sort();
547
548 // loop over all muons
549 iterator.Reset();
550 while((candidate = static_cast<Candidate*>(iterator.Next())))
551 {
552 const TLorentzVector &momentum = candidate->Momentum;
553 const TLorentzVector &position = candidate->Position;
554
555 pt = momentum.Pt();
556 cosTheta = TMath::Abs(momentum.CosTheta());
557 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
558 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
559 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
560
561 entry = static_cast<Muon*>(branch->NewEntry());
562
563 entry->SetBit(kIsReferenced);
564 entry->SetUniqueID(candidate->GetUniqueID());
565
566 entry->Eta = eta;
567 entry->Phi = momentum.Phi();
568 entry->PT = pt;
569
570 entry->T = position.T()*1.0E-3/c_light;
571
572 // Isolation variables
573
574 entry->IsolationVar = candidate->IsolationVar;
575 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr ;
576 entry->SumPtCharged = candidate->SumPtCharged ;
577 entry->SumPtNeutral = candidate->SumPtNeutral ;
578 entry->SumPtChargedPU = candidate->SumPtChargedPU ;
579 entry->SumPt = candidate->SumPt ;
580
581 entry->Charge = candidate->Charge;
582
583 entry->Particle = candidate->GetCandidates()->At(0);
584 }
585}
586
587//------------------------------------------------------------------------------
588
589void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
590{
591 TIter iterator(array);
592 Candidate *candidate = 0, *constituent = 0;
593 Jet *entry = 0;
594 Double_t pt, signPz, cosTheta, eta, rapidity;
595 Double_t ecalEnergy, hcalEnergy;
596 const Double_t c_light = 2.99792458E8;
597 Int_t i;
598
599 array->Sort();
600
601 // loop over all jets
602 iterator.Reset();
603 while((candidate = static_cast<Candidate*>(iterator.Next())))
604 {
605 TIter itConstituents(candidate->GetCandidates());
606
607 const TLorentzVector &momentum = candidate->Momentum;
608 const TLorentzVector &position = candidate->Position;
609
610 pt = momentum.Pt();
611 cosTheta = TMath::Abs(momentum.CosTheta());
612 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
613 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
614 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
615
616 entry = static_cast<Jet*>(branch->NewEntry());
617
618 entry->Eta = eta;
619 entry->Phi = momentum.Phi();
620 entry->PT = pt;
621
622 entry->T = position.T()*1.0E-3/c_light;
623
624 entry->Mass = momentum.M();
625
626 entry->Area = candidate->Area;
627
628 entry->DeltaEta = candidate->DeltaEta;
629 entry->DeltaPhi = candidate->DeltaPhi;
630
631 entry->Flavor = candidate->Flavor;
632 entry->FlavorAlgo = candidate->FlavorAlgo;
633 entry->FlavorPhys = candidate->FlavorPhys;
634
635 entry->BTag = candidate->BTag;
636
637 entry->BTagAlgo = candidate->BTagAlgo;
638 entry->BTagPhys = candidate->BTagPhys;
639
640 entry->TauTag = candidate->TauTag;
641
642 entry->Charge = candidate->Charge;
643
644 itConstituents.Reset();
645 entry->Constituents.Clear();
646 ecalEnergy = 0.0;
647 hcalEnergy = 0.0;
648 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
649 {
650 entry->Constituents.Add(constituent);
651 ecalEnergy += constituent->Eem;
652 hcalEnergy += constituent->Ehad;
653 }
654
655 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
656
657 //--- Pile-Up Jet ID variables ----
658
659 entry->NCharged = candidate->NCharged;
660 entry->NNeutrals = candidate->NNeutrals;
661 entry->Beta = candidate->Beta;
662 entry->BetaStar = candidate->BetaStar;
663 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
664 entry->PTD = candidate->PTD;
665
666 //--- Sub-structure variables ----
667
668 entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
669 entry->NSubJetsPruned = candidate->NSubJetsPruned;
670 entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
671
672 for(i = 0; i < 5; i++)
673 {
674 entry->FracPt[i] = candidate -> FracPt[i];
675 entry->Tau[i] = candidate -> Tau[i];
676 entry->TrimmedP4[i] = candidate -> TrimmedP4[i];
677 entry->PrunedP4[i] = candidate -> PrunedP4[i];
678 entry->SoftDroppedP4[i] = candidate -> SoftDroppedP4[i];
679 }
680
681 FillParticles(candidate, &entry->Particles);
682 }
683}
684
685//------------------------------------------------------------------------------
686
687void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
688{
689 Candidate *candidate = 0;
690 MissingET *entry = 0;
691
692 // get the first entry
693 if((candidate = static_cast<Candidate*>(array->At(0))))
694 {
695 const TLorentzVector &momentum = candidate->Momentum;
696
697 entry = static_cast<MissingET*>(branch->NewEntry());
698
699 entry->Eta = (-momentum).Eta();
700 entry->Phi = (-momentum).Phi();
701 entry->MET = momentum.Pt();
702 }
703}
704
705//------------------------------------------------------------------------------
706
707void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
708{
709 Candidate *candidate = 0;
710 ScalarHT *entry = 0;
711
712 // get the first entry
713 if((candidate = static_cast<Candidate*>(array->At(0))))
714 {
715 const TLorentzVector &momentum = candidate->Momentum;
716
717 entry = static_cast<ScalarHT*>(branch->NewEntry());
718
719 entry->HT = momentum.Pt();
720 }
721}
722
723//------------------------------------------------------------------------------
724
725void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
726{
727 TIter iterator(array);
728 Candidate *candidate = 0;
729 Rho *entry = 0;
730
731 // loop over all rho
732 iterator.Reset();
733 while((candidate = static_cast<Candidate*>(iterator.Next())))
734 {
735 const TLorentzVector &momentum = candidate->Momentum;
736
737 entry = static_cast<Rho*>(branch->NewEntry());
738
739 entry->Rho = momentum.E();
740 entry->Edges[0] = candidate->Edges[0];
741 entry->Edges[1] = candidate->Edges[1];
742 }
743}
744
745//------------------------------------------------------------------------------
746
747void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
748{
749 Candidate *candidate = 0;
750 Weight *entry = 0;
751
752 // get the first entry
753 if((candidate = static_cast<Candidate*>(array->At(0))))
754 {
755 const TLorentzVector &momentum = candidate->Momentum;
756
757 entry = static_cast<Weight*>(branch->NewEntry());
758
759 entry->Weight = momentum.E();
760 }
761}
762
763//------------------------------------------------------------------------------
764
765void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
766{
767 TIter iterator(array);
768 Candidate *candidate = 0;
769 HectorHit *entry = 0;
770
771 // loop over all roman pot hits
772 iterator.Reset();
773 while((candidate = static_cast<Candidate*>(iterator.Next())))
774 {
775 const TLorentzVector &position = candidate->Position;
776 const TLorentzVector &momentum = candidate->Momentum;
777
778 entry = static_cast<HectorHit*>(branch->NewEntry());
779
780 entry->E = momentum.E();
781
782 entry->Tx = momentum.Px();
783 entry->Ty = momentum.Py();
784
785 entry->T = position.T();
786
787 entry->X = position.X();
788 entry->Y = position.Y();
789 entry->S = position.Z();
790
791 entry->Particle = candidate->GetCandidates()->At(0);
792 }
793}
794
795//------------------------------------------------------------------------------
796
797void TreeWriter::Process()
798{
799 TBranchMap::iterator itBranchMap;
800 ExRootTreeBranch *branch;
801 TProcessMethod method;
802 TObjArray *array;
803
804 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
805 {
806 branch = itBranchMap->first;
807 method = itBranchMap->second.first;
808 array = itBranchMap->second.second;
809
810 (this->*method)(branch, array);
811 }
812}
813
814//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.