Fork me on GitHub

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

Last change on this file since 3a105e5 was 3a105e5, checked in by michele <michele.selvaggi@…>, 3 years ago

added first hit to track

  • Property mode set to 100644
File size: 29.5 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, m;
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 entry->Nclusters = candidate->Nclusters;
358 entry->dNdx = candidate->dNdx;
359
360 entry->ErrorP = candidate->ErrorP;
361 entry->ErrorPT = candidate->ErrorPT;
362
363 // diagonal covariance matrix terms
364 entry->ErrorD0 = candidate->ErrorD0;
365 entry->ErrorC = candidate->ErrorC;
366 entry->ErrorPhi = candidate->ErrorPhi;
367 entry->ErrorDZ = candidate->ErrorDZ;
368 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
369
370 // add some offdiagonal covariance matrix elements
371 entry->ErrorD0Phi = candidate->TrackCovariance(0,1)*1.e3;
372 entry->ErrorD0C = candidate->TrackCovariance(0,2);
373 entry->ErrorD0DZ = candidate->TrackCovariance(0,3)*1.e6;
374 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4)*1.e3;
375 entry->ErrorPhiC = candidate->TrackCovariance(1,2)*1.e-3;
376 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3)*1.e3;
377 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
378 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
379 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4)*1.e-3;
380 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4)*1.e3;
381
382 entry->Xd = candidate->Xd;
383 entry->Yd = candidate->Yd;
384 entry->Zd = candidate->Zd;
385
386 entry->XFirstHit = candidate->XFirstHit;
387 entry->YFirstHit = candidate->YFirstHit;
388 entry->ZFirstHit = candidate->ZFirstHit;
389
390 const TLorentzVector &momentum = candidate->Momentum;
391
392 pt = momentum.Pt();
393 p = momentum.P();
394 phi = momentum.Phi();
395 m = momentum.M();
396 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
397
398 cosTheta = TMath::Abs(momentum.CosTheta());
399 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
400 eta = (cosTheta == 1.0 ? signz * 999.9 : momentum.Eta());
401 rapidity = (cosTheta == 1.0 ? signz * 999.9 : momentum.Rapidity());
402
403 entry->P = p;
404 entry->PT = pt;
405 entry->Eta = eta;
406 entry->Phi = phi;
407 entry->CtgTheta = ctgTheta;
408 entry->C = candidate->C;
409 entry->Mass = m;
410
411 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
412 //const TLorentzVector &initialPosition = particle->Position;
413 const TLorentzVector &initialPosition = candidate->InitialPosition;
414
415 entry->X = initialPosition.X();
416 entry->Y = initialPosition.Y();
417 entry->Z = initialPosition.Z();
418 entry->T = initialPosition.T() * 1.0E-3 / c_light;
419 entry->ErrorT =candidate-> ErrorT * 1.0E-3 / c_light;
420
421 entry->Particle = particle;
422
423 entry->VertexIndex = candidate->ClusterIndex;
424 }
425}
426
427//------------------------------------------------------------------------------
428
429void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
430{
431 TIter iterator(array);
432 Candidate *candidate = 0;
433 Tower *entry = 0;
434 Double_t pt, signPz, cosTheta, eta, rapidity;
435 const Double_t c_light = 2.99792458E8;
436
437 // loop over all towers
438 iterator.Reset();
439 while((candidate = static_cast<Candidate *>(iterator.Next())))
440 {
441 const TLorentzVector &momentum = candidate->Momentum;
442 const TLorentzVector &position = candidate->Position;
443
444 pt = momentum.Pt();
445 cosTheta = TMath::Abs(momentum.CosTheta());
446 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
447 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
448 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
449
450 entry = static_cast<Tower *>(branch->NewEntry());
451
452 entry->SetBit(kIsReferenced);
453 entry->SetUniqueID(candidate->GetUniqueID());
454
455 entry->Eta = eta;
456 entry->Phi = momentum.Phi();
457 entry->ET = pt;
458 entry->E = momentum.E();
459 entry->Eem = candidate->Eem;
460 entry->Ehad = candidate->Ehad;
461 entry->Etrk = candidate->Etrk;
462 entry->Edges[0] = candidate->Edges[0];
463 entry->Edges[1] = candidate->Edges[1];
464 entry->Edges[2] = candidate->Edges[2];
465 entry->Edges[3] = candidate->Edges[3];
466
467 entry->T = position.T() * 1.0E-3 / c_light;
468 entry->NTimeHits = candidate->NTimeHits;
469
470 FillParticles(candidate, &entry->Particles);
471 }
472}
473
474//------------------------------------------------------------------------------
475
476void TreeWriter::ProcessParticleFlowCandidates(ExRootTreeBranch *branch, TObjArray *array)
477{
478
479 TIter iterator(array);
480 Candidate *candidate = 0;
481 Candidate *particle = 0;
482 ParticleFlowCandidate *entry = 0;
483 Double_t e, pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
484 const Double_t c_light = 2.99792458E8;
485
486 // loop over all tracks
487 iterator.Reset();
488 while((candidate = static_cast<Candidate *>(iterator.Next())))
489 {
490 const TLorentzVector &position = candidate->Position;
491
492 cosTheta = TMath::Abs(position.CosTheta());
493 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
494 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
495 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
496
497 entry = static_cast<ParticleFlowCandidate *>(branch->NewEntry());
498
499 entry->SetBit(kIsReferenced);
500 entry->SetUniqueID(candidate->GetUniqueID());
501
502 entry->PID = candidate->PID;
503
504 entry->Charge = candidate->Charge;
505
506 entry->EtaOuter = eta;
507 entry->PhiOuter = position.Phi();
508
509 entry->XOuter = position.X();
510 entry->YOuter = position.Y();
511 entry->ZOuter = position.Z();
512 entry->TOuter = position.T() * 1.0E-3 / c_light;
513
514 entry->L = candidate->L;
515
516 entry->D0 = candidate->D0;
517 entry->DZ = candidate->DZ;
518 entry->Nclusters = candidate->Nclusters;
519 entry->dNdx = candidate->dNdx;
520
521 entry->ErrorP = candidate->ErrorP;
522 entry->ErrorPT = candidate->ErrorPT;
523 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
524
525
526 // diagonal covariance matrix terms
527
528 entry->ErrorD0 = candidate->ErrorD0;
529 entry->ErrorC = candidate->ErrorC;
530 entry->ErrorPhi = candidate->ErrorPhi;
531 entry->ErrorDZ = candidate->ErrorDZ;
532 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
533
534 // add some offdiagonal covariance matrix elements
535 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
536 entry->ErrorD0C = candidate->TrackCovariance(0,2);
537 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
538 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
539 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
540 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
541 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
542 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
543 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
544 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
545
546 entry->Xd = candidate->Xd;
547 entry->Yd = candidate->Yd;
548 entry->Zd = candidate->Zd;
549
550 entry->XFirstHit = candidate->XFirstHit;
551 entry->YFirstHit = candidate->YFirstHit;
552 entry->ZFirstHit = candidate->ZFirstHit;
553
554 const TLorentzVector &momentum = candidate->Momentum;
555
556 e = momentum.E();
557 pt = momentum.Pt();
558 p = momentum.P();
559 phi = momentum.Phi();
560 m = momentum.M();
561 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
562
563 entry->E = e;
564 entry->P = p;
565 entry->PT = pt;
566 entry->Eta = eta;
567 entry->Phi = phi;
568 entry->CtgTheta = ctgTheta;
569 entry->C = candidate->C;
570 entry->Mass = m;
571
572 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
573 //const TLorentzVector &initialPosition = particle->Position;
574 const TLorentzVector &initialPosition = candidate->InitialPosition;
575
576 entry->X = initialPosition.X();
577 entry->Y = initialPosition.Y();
578 entry->Z = initialPosition.Z();
579 entry->T = initialPosition.T() * 1.0E-3 / c_light;
580 entry->ErrorT = candidate-> ErrorT * 1.0E-3 / c_light;
581
582 entry->VertexIndex = candidate->ClusterIndex;
583
584 entry->Eem = candidate->Eem;
585 entry->Ehad = candidate->Ehad;
586 entry->Etrk = candidate->Etrk;
587 entry->Edges[0] = candidate->Edges[0];
588 entry->Edges[1] = candidate->Edges[1];
589 entry->Edges[2] = candidate->Edges[2];
590 entry->Edges[3] = candidate->Edges[3];
591
592 //entry->T = position.T() * 1.0E-3 / c_light;
593 entry->NTimeHits = candidate->NTimeHits;
594
595 FillParticles(candidate, &entry->Particles);
596
597 }
598}
599
600//------------------------------------------------------------------------------
601
602void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
603{
604 TIter iterator(array);
605 Candidate *candidate = 0;
606 Photon *entry = 0;
607 Double_t pt, signPz, cosTheta, eta, rapidity;
608 const Double_t c_light = 2.99792458E8;
609
610 array->Sort();
611
612 // loop over all photons
613 iterator.Reset();
614 while((candidate = static_cast<Candidate *>(iterator.Next())))
615 {
616 TIter it1(candidate->GetCandidates());
617 const TLorentzVector &momentum = candidate->Momentum;
618 const TLorentzVector &position = candidate->Position;
619
620 pt = momentum.Pt();
621 cosTheta = TMath::Abs(momentum.CosTheta());
622 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
623 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
624 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
625
626 entry = static_cast<Photon *>(branch->NewEntry());
627
628 entry->Eta = eta;
629 entry->Phi = momentum.Phi();
630 entry->PT = pt;
631 entry->E = momentum.E();
632 entry->T = position.T() * 1.0E-3 / c_light;
633
634 // Isolation variables
635
636 entry->IsolationVar = candidate->IsolationVar;
637 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
638 entry->SumPtCharged = candidate->SumPtCharged;
639 entry->SumPtNeutral = candidate->SumPtNeutral;
640 entry->SumPtChargedPU = candidate->SumPtChargedPU;
641 entry->SumPt = candidate->SumPt;
642
643 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad / candidate->Eem : 999.9;
644
645 // 1: prompt -- 2: non prompt -- 3: fake
646 entry->Status = candidate->Status;
647
648 FillParticles(candidate, &entry->Particles);
649 }
650}
651
652//------------------------------------------------------------------------------
653
654void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
655{
656 TIter iterator(array);
657 Candidate *candidate = 0;
658 Electron *entry = 0;
659 Double_t pt, signPz, cosTheta, eta, rapidity;
660 const Double_t c_light = 2.99792458E8;
661
662 array->Sort();
663
664 // loop over all electrons
665 iterator.Reset();
666 while((candidate = static_cast<Candidate *>(iterator.Next())))
667 {
668 const TLorentzVector &momentum = candidate->Momentum;
669 const TLorentzVector &position = candidate->Position;
670
671 pt = momentum.Pt();
672 cosTheta = TMath::Abs(momentum.CosTheta());
673 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
674 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
675 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
676
677 entry = static_cast<Electron *>(branch->NewEntry());
678
679 entry->Eta = eta;
680 entry->Phi = momentum.Phi();
681 entry->PT = pt;
682
683 entry->T = position.T() * 1.0E-3 / c_light;
684
685 // displacement
686 entry->D0 = candidate->D0;
687 entry->ErrorD0 = candidate->ErrorD0;
688 entry->DZ = candidate->DZ;
689 entry->ErrorDZ = candidate->ErrorDZ;
690
691 // Isolation variables
692 entry->IsolationVar = candidate->IsolationVar;
693 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
694 entry->SumPtCharged = candidate->SumPtCharged;
695 entry->SumPtNeutral = candidate->SumPtNeutral;
696 entry->SumPtChargedPU = candidate->SumPtChargedPU;
697 entry->SumPt = candidate->SumPt;
698
699 entry->Charge = candidate->Charge;
700
701 entry->EhadOverEem = 0.0;
702
703 entry->Particle = candidate->GetCandidates()->At(0);
704 }
705}
706
707//------------------------------------------------------------------------------
708
709void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
710{
711 TIter iterator(array);
712 Candidate *candidate = 0;
713 Muon *entry = 0;
714 Double_t pt, signPz, cosTheta, eta, rapidity;
715
716 const Double_t c_light = 2.99792458E8;
717
718 array->Sort();
719
720 // loop over all muons
721 iterator.Reset();
722 while((candidate = static_cast<Candidate *>(iterator.Next())))
723 {
724 const TLorentzVector &momentum = candidate->Momentum;
725 const TLorentzVector &position = candidate->Position;
726
727 pt = momentum.Pt();
728 cosTheta = TMath::Abs(momentum.CosTheta());
729 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
730 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
731 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
732
733 entry = static_cast<Muon *>(branch->NewEntry());
734
735 entry->SetBit(kIsReferenced);
736 entry->SetUniqueID(candidate->GetUniqueID());
737
738 entry->Eta = eta;
739 entry->Phi = momentum.Phi();
740 entry->PT = pt;
741
742 entry->T = position.T() * 1.0E-3 / c_light;
743
744 // displacement
745 entry->D0 = candidate->D0;
746 entry->ErrorD0 = candidate->ErrorD0;
747 entry->DZ = candidate->DZ;
748 entry->ErrorDZ = candidate->ErrorDZ;
749
750 // Isolation variables
751
752 entry->IsolationVar = candidate->IsolationVar;
753 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
754 entry->SumPtCharged = candidate->SumPtCharged;
755 entry->SumPtNeutral = candidate->SumPtNeutral;
756 entry->SumPtChargedPU = candidate->SumPtChargedPU;
757 entry->SumPt = candidate->SumPt;
758
759 entry->Charge = candidate->Charge;
760
761 entry->Particle = candidate->GetCandidates()->At(0);
762 }
763}
764
765//------------------------------------------------------------------------------
766
767void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
768{
769 TIter iterator(array);
770 Candidate *candidate = 0, *constituent = 0;
771 Jet *entry = 0;
772 Double_t pt, signPz, cosTheta, eta, rapidity;
773 Double_t ecalEnergy, hcalEnergy;
774 const Double_t c_light = 2.99792458E8;
775 Int_t i;
776
777 array->Sort();
778
779 // loop over all jets
780 iterator.Reset();
781 while((candidate = static_cast<Candidate *>(iterator.Next())))
782 {
783 TIter itConstituents(candidate->GetCandidates());
784
785 const TLorentzVector &momentum = candidate->Momentum;
786 const TLorentzVector &position = candidate->Position;
787
788 pt = momentum.Pt();
789 cosTheta = TMath::Abs(momentum.CosTheta());
790 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
791 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
792 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
793
794 entry = static_cast<Jet *>(branch->NewEntry());
795
796 entry->Eta = eta;
797 entry->Phi = momentum.Phi();
798 entry->PT = pt;
799
800 entry->T = position.T() * 1.0E-3 / c_light;
801
802 entry->Mass = momentum.M();
803
804 entry->Area = candidate->Area;
805
806 entry->DeltaEta = candidate->DeltaEta;
807 entry->DeltaPhi = candidate->DeltaPhi;
808
809 entry->Flavor = candidate->Flavor;
810 entry->FlavorAlgo = candidate->FlavorAlgo;
811 entry->FlavorPhys = candidate->FlavorPhys;
812
813 entry->BTag = candidate->BTag;
814
815 entry->BTagAlgo = candidate->BTagAlgo;
816 entry->BTagPhys = candidate->BTagPhys;
817
818 entry->TauTag = candidate->TauTag;
819 entry->TauWeight = candidate->TauWeight;
820
821 entry->Charge = candidate->Charge;
822
823 itConstituents.Reset();
824 entry->Constituents.Clear();
825 ecalEnergy = 0.0;
826 hcalEnergy = 0.0;
827 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
828 {
829 entry->Constituents.Add(constituent);
830 ecalEnergy += constituent->Eem;
831 hcalEnergy += constituent->Ehad;
832 }
833
834 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy / ecalEnergy : 999.9;
835
836 //--- Pile-Up Jet ID variables ----
837
838 entry->NCharged = candidate->NCharged;
839 entry->NNeutrals = candidate->NNeutrals;
840
841 entry->NeutralEnergyFraction = candidate->NeutralEnergyFraction;
842 entry->ChargedEnergyFraction = candidate->ChargedEnergyFraction;
843 entry->Beta = candidate->Beta;
844 entry->BetaStar = candidate->BetaStar;
845 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
846 entry->PTD = candidate->PTD;
847
848 //--- Sub-structure variables ----
849
850 entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
851 entry->NSubJetsPruned = candidate->NSubJetsPruned;
852 entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
853
854 entry->SoftDroppedJet = candidate->SoftDroppedJet;
855 entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1;
856 entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
857
858 for(i = 0; i < 5; i++)
859 {
860 entry->FracPt[i] = candidate->FracPt[i];
861 entry->Tau[i] = candidate->Tau[i];
862 entry->TrimmedP4[i] = candidate->TrimmedP4[i];
863 entry->PrunedP4[i] = candidate->PrunedP4[i];
864 entry->SoftDroppedP4[i] = candidate->SoftDroppedP4[i];
865 }
866
867 //--- exclusive clustering variables ---
868 entry->ExclYmerge23 = candidate->ExclYmerge23;
869 entry->ExclYmerge34 = candidate->ExclYmerge34;
870 entry->ExclYmerge45 = candidate->ExclYmerge45;
871 entry->ExclYmerge56 = candidate->ExclYmerge56;
872
873 FillParticles(candidate, &entry->Particles);
874 }
875}
876
877//------------------------------------------------------------------------------
878
879void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
880{
881 Candidate *candidate = 0;
882 MissingET *entry = 0;
883
884 // get the first entry
885 if((candidate = static_cast<Candidate *>(array->At(0))))
886 {
887 const TLorentzVector &momentum = candidate->Momentum;
888
889 entry = static_cast<MissingET *>(branch->NewEntry());
890
891 entry->Eta = (-momentum).Eta();
892 entry->Phi = (-momentum).Phi();
893 entry->MET = momentum.Pt();
894 }
895}
896
897//------------------------------------------------------------------------------
898
899void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
900{
901 Candidate *candidate = 0;
902 ScalarHT *entry = 0;
903
904 // get the first entry
905 if((candidate = static_cast<Candidate *>(array->At(0))))
906 {
907 const TLorentzVector &momentum = candidate->Momentum;
908
909 entry = static_cast<ScalarHT *>(branch->NewEntry());
910
911 entry->HT = momentum.Pt();
912 }
913}
914
915//------------------------------------------------------------------------------
916
917void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
918{
919 TIter iterator(array);
920 Candidate *candidate = 0;
921 Rho *entry = 0;
922
923 // loop over all rho
924 iterator.Reset();
925 while((candidate = static_cast<Candidate *>(iterator.Next())))
926 {
927 const TLorentzVector &momentum = candidate->Momentum;
928
929 entry = static_cast<Rho *>(branch->NewEntry());
930
931 entry->Rho = momentum.E();
932 entry->Edges[0] = candidate->Edges[0];
933 entry->Edges[1] = candidate->Edges[1];
934 }
935}
936
937//------------------------------------------------------------------------------
938
939void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
940{
941 Candidate *candidate = 0;
942 Weight *entry = 0;
943
944 // get the first entry
945 if((candidate = static_cast<Candidate *>(array->At(0))))
946 {
947 const TLorentzVector &momentum = candidate->Momentum;
948
949 entry = static_cast<Weight *>(branch->NewEntry());
950
951 entry->Weight = momentum.E();
952 }
953}
954
955//------------------------------------------------------------------------------
956
957void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
958{
959 TIter iterator(array);
960 Candidate *candidate = 0;
961 HectorHit *entry = 0;
962
963 // loop over all roman pot hits
964 iterator.Reset();
965 while((candidate = static_cast<Candidate *>(iterator.Next())))
966 {
967 const TLorentzVector &position = candidate->Position;
968 const TLorentzVector &momentum = candidate->Momentum;
969
970 entry = static_cast<HectorHit *>(branch->NewEntry());
971
972 entry->E = momentum.E();
973
974 entry->Tx = momentum.Px();
975 entry->Ty = momentum.Py();
976
977 entry->T = position.T();
978
979 entry->X = position.X();
980 entry->Y = position.Y();
981 entry->S = position.Z();
982
983 entry->Particle = candidate->GetCandidates()->At(0);
984 }
985}
986
987//------------------------------------------------------------------------------
988
989void TreeWriter::Process()
990{
991 TBranchMap::iterator itBranchMap;
992 ExRootTreeBranch *branch;
993 TProcessMethod method;
994 TObjArray *array;
995
996 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
997 {
998 branch = itBranchMap->first;
999 method = itBranchMap->second.first;
1000 array = itBranchMap->second.second;
1001
1002 (this->*method)(branch, array);
1003 }
1004}
1005
1006//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.