Fork me on GitHub

source: git/modules/TreeWriter.cc@ 4265189

Timing
Last change on this file since 4265189 was 84a097e, checked in by Michele Selvaggi <michele.selvaggi@…>, 5 years ago

first iteration of adding timing in Delphes

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