Fork me on GitHub

source: git/modules/TreeWriter.cc@ 04290b1

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

add time error in Vertex class

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