Fork me on GitHub

source: git/modules/TreeWriter.cc@ 3e2bb2b

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

fixed vertex refs

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