Fork me on GitHub

source: git/modules/TreeWriter.cc@ 40e890c

Last change on this file since 40e890c was 17cd992, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

added track curvature inverse to Track and PFcandidate ouput

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