Fork me on GitHub

source: git/modules/TreeWriter.cc@ 0e6b14e

ImprovedOutputFile Timing
Last change on this file since 0e6b14e was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

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