Fork me on GitHub

source: git/modules/TreeWriter.cc@ 17d0ab8

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

Merge branch 'master' into photonId

Conflicts:

cards/CMS_PhaseII/CMS_PhaseII_200PU_v03.tcl

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