Fork me on GitHub

source: git/modules/TreeWriter.cc@ 2eaa9fd

Last change on this file since 2eaa9fd was dc883b4, checked in by Pavel Demin <pavel.demin@…>, 4 years ago

add info to TreeWriter

  • Property mode set to 100644
File size: 28.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/** \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[ParticleFlowCandidate::Class()] = &TreeWriter::ProcessParticleFlowCandidates;
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 param = GetParam("Info");
127 TString infoName;
128 Double_t infoValue;
129
130 size = param.GetSize();
131 for(i = 0; i < size / 2; ++i)
132 {
133 infoName = param[i * 2].GetString();
134 infoValue = param[i * 2 + 1].GetDouble();
135
136 AddInfo(infoName, infoValue);
137 }
138}
139
140//------------------------------------------------------------------------------
141
142void TreeWriter::Finish()
143{
144}
145
146//------------------------------------------------------------------------------
147
148void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
149{
150 TIter it1(candidate->GetCandidates());
151 it1.Reset();
152 array->Clear();
153
154 while((candidate = static_cast<Candidate *>(it1.Next())))
155 {
156 TIter it2(candidate->GetCandidates());
157
158 // particle
159 if(candidate->GetCandidates()->GetEntriesFast() == 0)
160 {
161 array->Add(candidate);
162 continue;
163 }
164
165 // track
166 candidate = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
167 if(candidate->GetCandidates()->GetEntriesFast() == 0)
168 {
169 array->Add(candidate);
170 continue;
171 }
172
173 // tower
174 it2.Reset();
175 while((candidate = static_cast<Candidate *>(it2.Next())))
176 {
177 array->Add(candidate->GetCandidates()->At(0));
178 }
179 }
180}
181
182//------------------------------------------------------------------------------
183
184void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
185{
186 TIter iterator(array);
187 Candidate *candidate = 0;
188 GenParticle *entry = 0;
189 Double_t pt, signPz, cosTheta, eta, rapidity;
190
191 const Double_t c_light = 2.99792458E8;
192
193 // loop over all particles
194 iterator.Reset();
195 while((candidate = static_cast<Candidate *>(iterator.Next())))
196 {
197 const TLorentzVector &momentum = candidate->Momentum;
198 const TLorentzVector &position = candidate->Position;
199
200 entry = static_cast<GenParticle *>(branch->NewEntry());
201
202 entry->SetBit(kIsReferenced);
203 entry->SetUniqueID(candidate->GetUniqueID());
204
205 pt = momentum.Pt();
206 cosTheta = TMath::Abs(momentum.CosTheta());
207 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
208 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
209 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
210
211 entry->PID = candidate->PID;
212
213 entry->Status = candidate->Status;
214 entry->IsPU = candidate->IsPU;
215
216 entry->M1 = candidate->M1;
217 entry->M2 = candidate->M2;
218
219 entry->D1 = candidate->D1;
220 entry->D2 = candidate->D2;
221
222 entry->Charge = candidate->Charge;
223 entry->Mass = candidate->Mass;
224
225 entry->E = momentum.E();
226 entry->Px = momentum.Px();
227 entry->Py = momentum.Py();
228 entry->Pz = momentum.Pz();
229
230 entry->Eta = eta;
231 entry->Phi = momentum.Phi();
232 entry->PT = pt;
233
234 entry->Rapidity = rapidity;
235
236 entry->X = position.X();
237 entry->Y = position.Y();
238 entry->Z = position.Z();
239 entry->T = position.T() * 1.0E-3 / c_light;
240 }
241}
242
243//------------------------------------------------------------------------------
244
245void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
246{
247 TIter iterator(array);
248 Candidate *candidate = 0, *constituent = 0;
249 Vertex *entry = 0;
250
251 const Double_t c_light = 2.99792458E8;
252
253 Double_t x, y, z, t, xError, yError, zError, tError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2;
254 UInt_t index, ndf;
255
256 CompBase *compare = Candidate::fgCompare;
257 Candidate::fgCompare = CompSumPT2<Candidate>::Instance();
258 array->Sort();
259 Candidate::fgCompare = compare;
260
261 // loop over all vertices
262 iterator.Reset();
263 while((candidate = static_cast<Candidate *>(iterator.Next())))
264 {
265
266 index = candidate->ClusterIndex;
267 ndf = candidate->ClusterNDF;
268 sigma = candidate->ClusterSigma;
269 sumPT2 = candidate->SumPT2;
270 btvSumPT2 = candidate->BTVSumPT2;
271 genDeltaZ = candidate->GenDeltaZ;
272 genSumPT2 = candidate->GenSumPT2;
273
274 x = candidate->Position.X();
275 y = candidate->Position.Y();
276 z = candidate->Position.Z();
277 t = candidate->Position.T() * 1.0E-3 / c_light;
278
279 xError = candidate->PositionError.X();
280 yError = candidate->PositionError.Y();
281 zError = candidate->PositionError.Z();
282 tError = candidate->PositionError.T() * 1.0E-3 / c_light;
283
284 entry = static_cast<Vertex *>(branch->NewEntry());
285
286 entry->Index = index;
287 entry->NDF = ndf;
288 entry->Sigma = sigma;
289 entry->SumPT2 = sumPT2;
290 entry->BTVSumPT2 = btvSumPT2;
291 entry->GenDeltaZ = genDeltaZ;
292 entry->GenSumPT2 = genSumPT2;
293
294 entry->X = x;
295 entry->Y = y;
296 entry->Z = z;
297 entry->T = t;
298
299 entry->ErrorX = xError;
300 entry->ErrorY = yError;
301 entry->ErrorZ = zError;
302 entry->ErrorT = tError;
303
304 TIter itConstituents(candidate->GetCandidates());
305 itConstituents.Reset();
306 entry->Constituents.Clear();
307 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
308 {
309 entry->Constituents.Add(constituent);
310 }
311 }
312}
313
314//------------------------------------------------------------------------------
315
316void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
317{
318 TIter iterator(array);
319 Candidate *candidate = 0;
320 Candidate *particle = 0;
321 Track *entry = 0;
322 Double_t pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi;
323 const Double_t c_light = 2.99792458E8;
324
325 // loop over all tracks
326 iterator.Reset();
327 while((candidate = static_cast<Candidate *>(iterator.Next())))
328 {
329 const TLorentzVector &position = candidate->Position;
330
331 cosTheta = TMath::Abs(position.CosTheta());
332 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
333 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
334 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
335
336 entry = static_cast<Track *>(branch->NewEntry());
337
338 entry->SetBit(kIsReferenced);
339 entry->SetUniqueID(candidate->GetUniqueID());
340
341 entry->PID = candidate->PID;
342
343 entry->Charge = candidate->Charge;
344
345 entry->EtaOuter = eta;
346 entry->PhiOuter = position.Phi();
347
348 entry->XOuter = position.X();
349 entry->YOuter = position.Y();
350 entry->ZOuter = position.Z();
351 entry->TOuter = position.T() * 1.0E-3 / c_light;
352
353 entry->L = candidate->L;
354
355 entry->D0 = candidate->D0;
356 entry->DZ = candidate->DZ;
357
358 entry->ErrorP = candidate->ErrorP;
359 entry->ErrorPT = candidate->ErrorPT;
360
361 // diagonal covariance matrix terms
362 entry->ErrorD0 = candidate->ErrorD0;
363 entry->ErrorC = candidate->ErrorC;
364 entry->ErrorPhi = candidate->ErrorPhi;
365 entry->ErrorDZ = candidate->ErrorDZ;
366 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
367
368 // add some offdiagonal covariance matrix elements
369 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
370 entry->ErrorD0C = candidate->TrackCovariance(0,2);
371 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
372 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
373 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
374 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
375 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
376 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
377 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
378 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
379
380 entry->Xd = candidate->Xd;
381 entry->Yd = candidate->Yd;
382 entry->Zd = candidate->Zd;
383
384 const TLorentzVector &momentum = candidate->Momentum;
385
386 pt = momentum.Pt();
387 p = momentum.P();
388 phi = momentum.Phi();
389 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
390
391 cosTheta = TMath::Abs(momentum.CosTheta());
392 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
393 eta = (cosTheta == 1.0 ? signz * 999.9 : momentum.Eta());
394 rapidity = (cosTheta == 1.0 ? signz * 999.9 : momentum.Rapidity());
395
396 entry->P = p;
397 entry->PT = pt;
398 entry->Eta = eta;
399 entry->Phi = phi;
400 entry->CtgTheta = ctgTheta;
401 entry->C = candidate->C;
402
403 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
404 const TLorentzVector &initialPosition = particle->Position;
405
406 entry->X = initialPosition.X();
407 entry->Y = initialPosition.Y();
408 entry->Z = initialPosition.Z();
409 entry->T = initialPosition.T() * 1.0E-3 / c_light;
410
411 entry->Particle = particle;
412
413 entry->VertexIndex = candidate->ClusterIndex;
414 }
415}
416
417//------------------------------------------------------------------------------
418
419void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
420{
421 TIter iterator(array);
422 Candidate *candidate = 0;
423 Tower *entry = 0;
424 Double_t pt, signPz, cosTheta, eta, rapidity;
425 const Double_t c_light = 2.99792458E8;
426
427 // loop over all towers
428 iterator.Reset();
429 while((candidate = static_cast<Candidate *>(iterator.Next())))
430 {
431 const TLorentzVector &momentum = candidate->Momentum;
432 const TLorentzVector &position = candidate->Position;
433
434 pt = momentum.Pt();
435 cosTheta = TMath::Abs(momentum.CosTheta());
436 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
437 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
438 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
439
440 entry = static_cast<Tower *>(branch->NewEntry());
441
442 entry->SetBit(kIsReferenced);
443 entry->SetUniqueID(candidate->GetUniqueID());
444
445 entry->Eta = eta;
446 entry->Phi = momentum.Phi();
447 entry->ET = pt;
448 entry->E = momentum.E();
449 entry->Eem = candidate->Eem;
450 entry->Ehad = candidate->Ehad;
451 entry->Edges[0] = candidate->Edges[0];
452 entry->Edges[1] = candidate->Edges[1];
453 entry->Edges[2] = candidate->Edges[2];
454 entry->Edges[3] = candidate->Edges[3];
455
456 entry->T = position.T() * 1.0E-3 / c_light;
457 entry->NTimeHits = candidate->NTimeHits;
458
459 FillParticles(candidate, &entry->Particles);
460 }
461}
462
463//------------------------------------------------------------------------------
464
465void TreeWriter::ProcessParticleFlowCandidates(ExRootTreeBranch *branch, TObjArray *array)
466{
467
468 TIter iterator(array);
469 Candidate *candidate = 0;
470 Candidate *particle = 0;
471 ParticleFlowCandidate *entry = 0;
472 Double_t e, pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi;
473 const Double_t c_light = 2.99792458E8;
474
475 // loop over all tracks
476 iterator.Reset();
477 while((candidate = static_cast<Candidate *>(iterator.Next())))
478 {
479 const TLorentzVector &position = candidate->Position;
480
481 cosTheta = TMath::Abs(position.CosTheta());
482 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
483 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
484 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
485
486 entry = static_cast<ParticleFlowCandidate *>(branch->NewEntry());
487
488 entry->SetBit(kIsReferenced);
489 entry->SetUniqueID(candidate->GetUniqueID());
490
491 entry->PID = candidate->PID;
492
493 entry->Charge = candidate->Charge;
494
495 entry->EtaOuter = eta;
496 entry->PhiOuter = position.Phi();
497
498 entry->XOuter = position.X();
499 entry->YOuter = position.Y();
500 entry->ZOuter = position.Z();
501 entry->TOuter = position.T() * 1.0E-3 / c_light;
502
503 entry->L = candidate->L;
504
505 entry->D0 = candidate->D0;
506 entry->DZ = candidate->DZ;
507
508 entry->ErrorP = candidate->ErrorP;
509 entry->ErrorPT = candidate->ErrorPT;
510 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
511
512
513 // diagonal covariance matrix terms
514
515 entry->ErrorD0 = candidate->ErrorD0;
516 entry->ErrorC = candidate->ErrorC;
517 entry->ErrorPhi = candidate->ErrorPhi;
518 entry->ErrorDZ = candidate->ErrorDZ;
519 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
520
521 // add some offdiagonal covariance matrix elements
522 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
523 entry->ErrorD0C = candidate->TrackCovariance(0,2);
524 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
525 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
526 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
527 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
528 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
529 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
530 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
531 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
532
533 entry->Xd = candidate->Xd;
534 entry->Yd = candidate->Yd;
535 entry->Zd = candidate->Zd;
536
537 const TLorentzVector &momentum = candidate->Momentum;
538
539 e = momentum.E();
540 pt = momentum.Pt();
541 p = momentum.P();
542 phi = momentum.Phi();
543 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
544
545 entry->E = e;
546 entry->P = p;
547 entry->PT = pt;
548 entry->Eta = eta;
549 entry->Phi = phi;
550 entry->CtgTheta = ctgTheta;
551 entry->C = candidate->C;
552
553 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
554 const TLorentzVector &initialPosition = particle->Position;
555
556 entry->X = initialPosition.X();
557 entry->Y = initialPosition.Y();
558 entry->Z = initialPosition.Z();
559 entry->T = initialPosition.T() * 1.0E-3 / c_light;
560
561 entry->VertexIndex = candidate->ClusterIndex;
562
563 entry->Eem = candidate->Eem;
564 entry->Ehad = candidate->Ehad;
565 entry->Edges[0] = candidate->Edges[0];
566 entry->Edges[1] = candidate->Edges[1];
567 entry->Edges[2] = candidate->Edges[2];
568 entry->Edges[3] = candidate->Edges[3];
569
570 entry->T = position.T() * 1.0E-3 / c_light;
571 entry->NTimeHits = candidate->NTimeHits;
572
573 FillParticles(candidate, &entry->Particles);
574
575 }
576}
577
578//------------------------------------------------------------------------------
579
580void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
581{
582 TIter iterator(array);
583 Candidate *candidate = 0;
584 Photon *entry = 0;
585 Double_t pt, signPz, cosTheta, eta, rapidity;
586 const Double_t c_light = 2.99792458E8;
587
588 array->Sort();
589
590 // loop over all photons
591 iterator.Reset();
592 while((candidate = static_cast<Candidate *>(iterator.Next())))
593 {
594 TIter it1(candidate->GetCandidates());
595 const TLorentzVector &momentum = candidate->Momentum;
596 const TLorentzVector &position = candidate->Position;
597
598 pt = momentum.Pt();
599 cosTheta = TMath::Abs(momentum.CosTheta());
600 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
601 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
602 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
603
604 entry = static_cast<Photon *>(branch->NewEntry());
605
606 entry->Eta = eta;
607 entry->Phi = momentum.Phi();
608 entry->PT = pt;
609 entry->E = momentum.E();
610 entry->T = position.T() * 1.0E-3 / c_light;
611
612 // Isolation variables
613
614 entry->IsolationVar = candidate->IsolationVar;
615 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
616 entry->SumPtCharged = candidate->SumPtCharged;
617 entry->SumPtNeutral = candidate->SumPtNeutral;
618 entry->SumPtChargedPU = candidate->SumPtChargedPU;
619 entry->SumPt = candidate->SumPt;
620
621 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad / candidate->Eem : 999.9;
622
623 // 1: prompt -- 2: non prompt -- 3: fake
624 entry->Status = candidate->Status;
625
626 FillParticles(candidate, &entry->Particles);
627 }
628}
629
630//------------------------------------------------------------------------------
631
632void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
633{
634 TIter iterator(array);
635 Candidate *candidate = 0;
636 Electron *entry = 0;
637 Double_t pt, signPz, cosTheta, eta, rapidity;
638 const Double_t c_light = 2.99792458E8;
639
640 array->Sort();
641
642 // loop over all electrons
643 iterator.Reset();
644 while((candidate = static_cast<Candidate *>(iterator.Next())))
645 {
646 const TLorentzVector &momentum = candidate->Momentum;
647 const TLorentzVector &position = candidate->Position;
648
649 pt = momentum.Pt();
650 cosTheta = TMath::Abs(momentum.CosTheta());
651 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
652 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
653 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
654
655 entry = static_cast<Electron *>(branch->NewEntry());
656
657 entry->Eta = eta;
658 entry->Phi = momentum.Phi();
659 entry->PT = pt;
660
661 entry->T = position.T() * 1.0E-3 / c_light;
662
663 // displacement
664 entry->D0 = candidate->D0;
665 entry->ErrorD0 = candidate->ErrorD0;
666 entry->DZ = candidate->DZ;
667 entry->ErrorDZ = candidate->ErrorDZ;
668
669 // Isolation variables
670 entry->IsolationVar = candidate->IsolationVar;
671 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
672 entry->SumPtCharged = candidate->SumPtCharged;
673 entry->SumPtNeutral = candidate->SumPtNeutral;
674 entry->SumPtChargedPU = candidate->SumPtChargedPU;
675 entry->SumPt = candidate->SumPt;
676
677 entry->Charge = candidate->Charge;
678
679 entry->EhadOverEem = 0.0;
680
681 entry->Particle = candidate->GetCandidates()->At(0);
682 }
683}
684
685//------------------------------------------------------------------------------
686
687void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
688{
689 TIter iterator(array);
690 Candidate *candidate = 0;
691 Muon *entry = 0;
692 Double_t pt, signPz, cosTheta, eta, rapidity;
693
694 const Double_t c_light = 2.99792458E8;
695
696 array->Sort();
697
698 // loop over all muons
699 iterator.Reset();
700 while((candidate = static_cast<Candidate *>(iterator.Next())))
701 {
702 const TLorentzVector &momentum = candidate->Momentum;
703 const TLorentzVector &position = candidate->Position;
704
705 pt = momentum.Pt();
706 cosTheta = TMath::Abs(momentum.CosTheta());
707 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
708 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
709 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
710
711 entry = static_cast<Muon *>(branch->NewEntry());
712
713 entry->SetBit(kIsReferenced);
714 entry->SetUniqueID(candidate->GetUniqueID());
715
716 entry->Eta = eta;
717 entry->Phi = momentum.Phi();
718 entry->PT = pt;
719
720 entry->T = position.T() * 1.0E-3 / c_light;
721
722 // displacement
723 entry->D0 = candidate->D0;
724 entry->ErrorD0 = candidate->ErrorD0;
725 entry->DZ = candidate->DZ;
726 entry->ErrorDZ = candidate->ErrorDZ;
727
728 // Isolation variables
729
730 entry->IsolationVar = candidate->IsolationVar;
731 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
732 entry->SumPtCharged = candidate->SumPtCharged;
733 entry->SumPtNeutral = candidate->SumPtNeutral;
734 entry->SumPtChargedPU = candidate->SumPtChargedPU;
735 entry->SumPt = candidate->SumPt;
736
737 entry->Charge = candidate->Charge;
738
739 entry->Particle = candidate->GetCandidates()->At(0);
740 }
741}
742
743//------------------------------------------------------------------------------
744
745void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
746{
747 TIter iterator(array);
748 Candidate *candidate = 0, *constituent = 0;
749 Jet *entry = 0;
750 Double_t pt, signPz, cosTheta, eta, rapidity;
751 Double_t ecalEnergy, hcalEnergy;
752 const Double_t c_light = 2.99792458E8;
753 Int_t i;
754
755 array->Sort();
756
757 // loop over all jets
758 iterator.Reset();
759 while((candidate = static_cast<Candidate *>(iterator.Next())))
760 {
761 TIter itConstituents(candidate->GetCandidates());
762
763 const TLorentzVector &momentum = candidate->Momentum;
764 const TLorentzVector &position = candidate->Position;
765
766 pt = momentum.Pt();
767 cosTheta = TMath::Abs(momentum.CosTheta());
768 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
769 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
770 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
771
772 entry = static_cast<Jet *>(branch->NewEntry());
773
774 entry->Eta = eta;
775 entry->Phi = momentum.Phi();
776 entry->PT = pt;
777
778 entry->T = position.T() * 1.0E-3 / c_light;
779
780 entry->Mass = momentum.M();
781
782 entry->Area = candidate->Area;
783
784 entry->DeltaEta = candidate->DeltaEta;
785 entry->DeltaPhi = candidate->DeltaPhi;
786
787 entry->Flavor = candidate->Flavor;
788 entry->FlavorAlgo = candidate->FlavorAlgo;
789 entry->FlavorPhys = candidate->FlavorPhys;
790
791 entry->BTag = candidate->BTag;
792
793 entry->BTagAlgo = candidate->BTagAlgo;
794 entry->BTagPhys = candidate->BTagPhys;
795
796 entry->TauTag = candidate->TauTag;
797 entry->TauWeight = candidate->TauWeight;
798
799 entry->Charge = candidate->Charge;
800
801 itConstituents.Reset();
802 entry->Constituents.Clear();
803 ecalEnergy = 0.0;
804 hcalEnergy = 0.0;
805 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
806 {
807 entry->Constituents.Add(constituent);
808 ecalEnergy += constituent->Eem;
809 hcalEnergy += constituent->Ehad;
810 }
811
812 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy / ecalEnergy : 999.9;
813
814 //--- Pile-Up Jet ID variables ----
815
816 entry->NCharged = candidate->NCharged;
817 entry->NNeutrals = candidate->NNeutrals;
818
819 entry->NeutralEnergyFraction = candidate->NeutralEnergyFraction;
820 entry->ChargedEnergyFraction = candidate->ChargedEnergyFraction;
821 entry->Beta = candidate->Beta;
822 entry->BetaStar = candidate->BetaStar;
823 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
824 entry->PTD = candidate->PTD;
825
826 //--- Sub-structure variables ----
827
828 entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
829 entry->NSubJetsPruned = candidate->NSubJetsPruned;
830 entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
831
832 entry->SoftDroppedJet = candidate->SoftDroppedJet;
833 entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1;
834 entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
835
836 for(i = 0; i < 5; i++)
837 {
838 entry->FracPt[i] = candidate->FracPt[i];
839 entry->Tau[i] = candidate->Tau[i];
840 entry->TrimmedP4[i] = candidate->TrimmedP4[i];
841 entry->PrunedP4[i] = candidate->PrunedP4[i];
842 entry->SoftDroppedP4[i] = candidate->SoftDroppedP4[i];
843 }
844
845 //--- exclusive clustering variables ---
846 entry->ExclYmerge23 = candidate->ExclYmerge23;
847 entry->ExclYmerge34 = candidate->ExclYmerge34;
848 entry->ExclYmerge45 = candidate->ExclYmerge45;
849 entry->ExclYmerge56 = candidate->ExclYmerge56;
850
851 FillParticles(candidate, &entry->Particles);
852 }
853}
854
855//------------------------------------------------------------------------------
856
857void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
858{
859 Candidate *candidate = 0;
860 MissingET *entry = 0;
861
862 // get the first entry
863 if((candidate = static_cast<Candidate *>(array->At(0))))
864 {
865 const TLorentzVector &momentum = candidate->Momentum;
866
867 entry = static_cast<MissingET *>(branch->NewEntry());
868
869 entry->Eta = (-momentum).Eta();
870 entry->Phi = (-momentum).Phi();
871 entry->MET = momentum.Pt();
872 }
873}
874
875//------------------------------------------------------------------------------
876
877void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
878{
879 Candidate *candidate = 0;
880 ScalarHT *entry = 0;
881
882 // get the first entry
883 if((candidate = static_cast<Candidate *>(array->At(0))))
884 {
885 const TLorentzVector &momentum = candidate->Momentum;
886
887 entry = static_cast<ScalarHT *>(branch->NewEntry());
888
889 entry->HT = momentum.Pt();
890 }
891}
892
893//------------------------------------------------------------------------------
894
895void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
896{
897 TIter iterator(array);
898 Candidate *candidate = 0;
899 Rho *entry = 0;
900
901 // loop over all rho
902 iterator.Reset();
903 while((candidate = static_cast<Candidate *>(iterator.Next())))
904 {
905 const TLorentzVector &momentum = candidate->Momentum;
906
907 entry = static_cast<Rho *>(branch->NewEntry());
908
909 entry->Rho = momentum.E();
910 entry->Edges[0] = candidate->Edges[0];
911 entry->Edges[1] = candidate->Edges[1];
912 }
913}
914
915//------------------------------------------------------------------------------
916
917void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
918{
919 Candidate *candidate = 0;
920 Weight *entry = 0;
921
922 // get the first entry
923 if((candidate = static_cast<Candidate *>(array->At(0))))
924 {
925 const TLorentzVector &momentum = candidate->Momentum;
926
927 entry = static_cast<Weight *>(branch->NewEntry());
928
929 entry->Weight = momentum.E();
930 }
931}
932
933//------------------------------------------------------------------------------
934
935void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
936{
937 TIter iterator(array);
938 Candidate *candidate = 0;
939 HectorHit *entry = 0;
940
941 // loop over all roman pot hits
942 iterator.Reset();
943 while((candidate = static_cast<Candidate *>(iterator.Next())))
944 {
945 const TLorentzVector &position = candidate->Position;
946 const TLorentzVector &momentum = candidate->Momentum;
947
948 entry = static_cast<HectorHit *>(branch->NewEntry());
949
950 entry->E = momentum.E();
951
952 entry->Tx = momentum.Px();
953 entry->Ty = momentum.Py();
954
955 entry->T = position.T();
956
957 entry->X = position.X();
958 entry->Y = position.Y();
959 entry->S = position.Z();
960
961 entry->Particle = candidate->GetCandidates()->At(0);
962 }
963}
964
965//------------------------------------------------------------------------------
966
967void TreeWriter::Process()
968{
969 TBranchMap::iterator itBranchMap;
970 ExRootTreeBranch *branch;
971 TProcessMethod method;
972 TObjArray *array;
973
974 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
975 {
976 branch = itBranchMap->first;
977 method = itBranchMap->second.first;
978 array = itBranchMap->second.second;
979
980 (this->*method)(branch, array);
981 }
982}
983
984//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.