Fork me on GitHub

source: git/modules/TreeWriter.cc@ 4e8e72b

Last change on this file since 4e8e72b was 4e8e72b, checked in by Pavel Demin <pavel-demin@…>, 3 years ago

deduplicate particles in TreeWriter::FillParticles

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