Fork me on GitHub

source: git/modules/TreeWriter.cc@ 8f7db23

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 8f7db23 was 8f7db23, checked in by pavel <pavel@…>, 11 years ago

add HectorHit class and corresponding ROOT tree branch

  • Property mode set to 100644
File size: 18.3 KB
Line 
1
2/** \class TreeWriter
3 *
4 * Fills ROOT tree branches.
5 *
6 * $Date$
7 * $Revision$
8 *
9 *
10 * \author P. Demin - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include "modules/TreeWriter.h"
15
16#include "classes/DelphesClasses.h"
17#include "classes/DelphesFactory.h"
18#include "classes/DelphesFormula.h"
19
20#include "ExRootAnalysis/ExRootResult.h"
21#include "ExRootAnalysis/ExRootFilter.h"
22#include "ExRootAnalysis/ExRootClassifier.h"
23#include "ExRootAnalysis/ExRootTreeBranch.h"
24
25#include "TROOT.h"
26#include "TMath.h"
27#include "TString.h"
28#include "TFormula.h"
29#include "TRandom3.h"
30#include "TObjArray.h"
31#include "TDatabasePDG.h"
32#include "TLorentzVector.h"
33
34#include <algorithm>
35#include <stdexcept>
36#include <iostream>
37#include <sstream>
38
39using namespace std;
40
41//------------------------------------------------------------------------------
42
43TreeWriter::TreeWriter()
44{
45}
46
47//------------------------------------------------------------------------------
48
49TreeWriter::~TreeWriter()
50{
51}
52
53//------------------------------------------------------------------------------
54
55void TreeWriter::Init()
56{
57 fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
58 fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
59 fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
60 fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
61 fClassMap[Photon::Class()] = &TreeWriter::ProcessPhotons;
62 fClassMap[Electron::Class()] = &TreeWriter::ProcessElectrons;
63 fClassMap[Muon::Class()] = &TreeWriter::ProcessMuons;
64 fClassMap[Jet::Class()] = &TreeWriter::ProcessJets;
65 fClassMap[MissingET::Class()] = &TreeWriter::ProcessMissingET;
66 fClassMap[ScalarHT::Class()] = &TreeWriter::ProcessScalarHT;
67 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
68 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
69 fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
70
71 TBranchMap::iterator itBranchMap;
72 map< TClass *, TProcessMethod >::iterator itClassMap;
73
74 // read branch configuration and
75 // import array with output from filter/classifier/jetfinder modules
76
77 ExRootConfParam param = GetParam("Branch");
78 Long_t i, size;
79 TString branchName, branchClassName, branchInputArray;
80 TClass *branchClass;
81 TObjArray *array;
82 ExRootTreeBranch *branch;
83
84 size = param.GetSize();
85 for(i = 0; i < size/3; ++i)
86 {
87 branchInputArray = param[i*3].GetString();
88 branchName = param[i*3 + 1].GetString();
89 branchClassName = param[i*3 + 2].GetString();
90
91 branchClass = gROOT->GetClass(branchClassName);
92
93 if(!branchClass)
94 {
95 cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
96 continue;
97 }
98
99 itClassMap = fClassMap.find(branchClass);
100 if(itClassMap == fClassMap.end())
101 {
102 cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
103 continue;
104 }
105
106 array = ImportArray(branchInputArray);
107 branch = NewBranch(branchName, branchClass);
108
109 fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
110 }
111
112}
113
114//------------------------------------------------------------------------------
115
116void TreeWriter::Finish()
117{
118}
119
120//------------------------------------------------------------------------------
121
122void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
123{
124 TIter it1(candidate->GetCandidates());
125 it1.Reset();
126 array->Clear();
127 while((candidate = static_cast<Candidate*>(it1.Next())))
128 {
129 TIter it2(candidate->GetCandidates());
130
131 // particle
132 if(candidate->GetCandidates()->GetEntriesFast() == 0)
133 {
134 array->Add(candidate);
135 continue;
136 }
137
138 // track
139 candidate = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
140 if(candidate->GetCandidates()->GetEntriesFast() == 0)
141 {
142 array->Add(candidate);
143 continue;
144 }
145
146 // tower
147 it2.Reset();
148 while((candidate = static_cast<Candidate*>(it2.Next())))
149 {
150 array->Add(candidate->GetCandidates()->At(0));
151 }
152 }
153}
154
155//------------------------------------------------------------------------------
156
157void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
158{
159 TIter iterator(array);
160 Candidate *candidate = 0;
161 GenParticle *entry = 0;
162 Double_t pt, signPz, cosTheta, eta, rapidity;
163
164 const Double_t c_light = 2.99792458E8;
165
166 // loop over all particles
167 iterator.Reset();
168 while((candidate = static_cast<Candidate*>(iterator.Next())))
169 {
170 const TLorentzVector &momentum = candidate->Momentum;
171 const TLorentzVector &position = candidate->Position;
172
173 entry = static_cast<GenParticle*>(branch->NewEntry());
174
175 entry->SetBit(kIsReferenced);
176 entry->SetUniqueID(candidate->GetUniqueID());
177
178 pt = momentum.Pt();
179 cosTheta = TMath::Abs(momentum.CosTheta());
180 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
181 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
182 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
183
184 entry->PID = candidate->PID;
185
186 entry->Status = candidate->Status;
187 entry->IsPU = candidate->IsPU;
188
189 entry->M1 = candidate->M1;
190 entry->M2 = candidate->M2;
191
192 entry->D1 = candidate->D1;
193 entry->D2 = candidate->D2;
194
195 entry->Charge = candidate->Charge;
196 entry->Mass = candidate->Mass;
197
198 entry->E = momentum.E();
199 entry->Px = momentum.Px();
200 entry->Py = momentum.Py();
201 entry->Pz = momentum.Pz();
202
203 entry->Eta = eta;
204 entry->Phi = momentum.Phi();
205 entry->PT = pt;
206
207 entry->Rapidity = rapidity;
208
209 entry->X = position.X();
210 entry->Y = position.Y();
211 entry->Z = position.Z();
212 entry->T = position.T()*1.0E-3/c_light;
213 }
214}
215
216//------------------------------------------------------------------------------
217
218void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
219{
220 TIter iterator(array);
221 Candidate *candidate = 0;
222 Vertex *entry = 0;
223
224 const Double_t c_light = 2.99792458E8;
225
226 // loop over all vertices
227 iterator.Reset();
228 while((candidate = static_cast<Candidate*>(iterator.Next())))
229 {
230 const TLorentzVector &position = candidate->Position;
231
232 entry = static_cast<Vertex*>(branch->NewEntry());
233
234 entry->X = position.X();
235 entry->Y = position.Y();
236 entry->Z = position.Z();
237 entry->T = position.T()*1.0E-3/c_light;
238 }
239}
240
241//------------------------------------------------------------------------------
242
243void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
244{
245 TIter iterator(array);
246 Candidate *candidate = 0;
247 Candidate *particle = 0;
248 Track *entry = 0;
249 Double_t pt, signz, cosTheta, eta, rapidity;
250 const Double_t c_light = 2.99792458E8;
251
252 // loop over all tracks
253 iterator.Reset();
254 while((candidate = static_cast<Candidate*>(iterator.Next())))
255 {
256 const TLorentzVector &position = candidate->Position;
257
258 cosTheta = TMath::Abs(position.CosTheta());
259 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
260 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
261 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
262
263 entry = static_cast<Track*>(branch->NewEntry());
264
265 entry->SetBit(kIsReferenced);
266 entry->SetUniqueID(candidate->GetUniqueID());
267
268 entry->PID = candidate->PID;
269
270 entry->Charge = candidate->Charge;
271
272 entry->EtaOuter = eta;
273 entry->PhiOuter = position.Phi();
274
275 entry->XOuter = position.X();
276 entry->YOuter = position.Y();
277 entry->ZOuter = position.Z();
278 entry->TOuter = position.T()*1.0E-3/c_light;
279
280 const TLorentzVector &momentum = candidate->Momentum;
281
282 pt = momentum.Pt();
283 cosTheta = TMath::Abs(momentum.CosTheta());
284 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
285 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
286 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
287
288 entry->Eta = eta;
289 entry->Phi = momentum.Phi();
290 entry->PT = pt;
291
292 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
293 const TLorentzVector &initialPosition = particle->Position;
294
295 entry->X = initialPosition.X();
296 entry->Y = initialPosition.Y();
297 entry->Z = initialPosition.Z();
298 entry->T = initialPosition.T()*1.0E-3/c_light;
299
300 entry->Particle = particle;
301 }
302}
303
304//------------------------------------------------------------------------------
305
306void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
307{
308 TIter iterator(array);
309 Candidate *candidate = 0;
310 Tower *entry = 0;
311 Double_t pt, signPz, cosTheta, eta, rapidity;
312 const Double_t c_light = 2.99792458E8;
313
314 // loop over all towers
315 iterator.Reset();
316 while((candidate = static_cast<Candidate*>(iterator.Next())))
317 {
318 const TLorentzVector &momentum = candidate->Momentum;
319 const TLorentzVector &position = candidate->Position;
320
321 pt = momentum.Pt();
322 cosTheta = TMath::Abs(momentum.CosTheta());
323 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
324 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
325 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
326
327 entry = static_cast<Tower*>(branch->NewEntry());
328
329 entry->SetBit(kIsReferenced);
330 entry->SetUniqueID(candidate->GetUniqueID());
331
332 entry->Eta = eta;
333 entry->Phi = momentum.Phi();
334 entry->ET = pt;
335 entry->E = momentum.E();
336 entry->Eem = candidate->Eem;
337 entry->Ehad = candidate->Ehad;
338 entry->Edges[0] = candidate->Edges[0];
339 entry->Edges[1] = candidate->Edges[1];
340 entry->Edges[2] = candidate->Edges[2];
341 entry->Edges[3] = candidate->Edges[3];
342
343 entry->T = position.T()*1.0E-3/c_light;
344
345 FillParticles(candidate, &entry->Particles);
346 }
347}
348
349//------------------------------------------------------------------------------
350
351void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
352{
353 TIter iterator(array);
354 Candidate *candidate = 0;
355 Photon *entry = 0;
356 Double_t pt, signPz, cosTheta, eta, rapidity;
357 const Double_t c_light = 2.99792458E8;
358
359 array->Sort();
360
361 // loop over all photons
362 iterator.Reset();
363 while((candidate = static_cast<Candidate*>(iterator.Next())))
364 {
365 TIter it1(candidate->GetCandidates());
366 const TLorentzVector &momentum = candidate->Momentum;
367 const TLorentzVector &position = candidate->Position;
368
369
370 pt = momentum.Pt();
371 cosTheta = TMath::Abs(momentum.CosTheta());
372 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
373 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
374 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
375
376 entry = static_cast<Photon*>(branch->NewEntry());
377
378 entry->Eta = eta;
379 entry->Phi = momentum.Phi();
380 entry->PT = pt;
381 entry->E = momentum.E();
382
383 entry->T = position.T()*1.0E-3/c_light;
384
385 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
386
387 FillParticles(candidate, &entry->Particles);
388 }
389}
390
391//------------------------------------------------------------------------------
392
393void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
394{
395 TIter iterator(array);
396 Candidate *candidate = 0;
397 Electron *entry = 0;
398 Double_t pt, signPz, cosTheta, eta, rapidity;
399 const Double_t c_light = 2.99792458E8;
400
401 array->Sort();
402
403 // loop over all electrons
404 iterator.Reset();
405 while((candidate = static_cast<Candidate*>(iterator.Next())))
406 {
407 const TLorentzVector &momentum = candidate->Momentum;
408 const TLorentzVector &position = candidate->Position;
409
410 pt = momentum.Pt();
411 cosTheta = TMath::Abs(momentum.CosTheta());
412 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
413 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
414 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
415
416 entry = static_cast<Electron*>(branch->NewEntry());
417
418 entry->Eta = eta;
419 entry->Phi = momentum.Phi();
420 entry->PT = pt;
421
422 entry->T = position.T()*1.0E-3/c_light;
423
424 entry->Charge = candidate->Charge;
425
426 entry->EhadOverEem = 0.0;
427
428 entry->Particle = candidate->GetCandidates()->At(0);
429 }
430}
431
432//------------------------------------------------------------------------------
433
434void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
435{
436 TIter iterator(array);
437 Candidate *candidate = 0;
438 Muon *entry = 0;
439 Double_t pt, signPz, cosTheta, eta, rapidity;
440
441 const Double_t c_light = 2.99792458E8;
442
443 array->Sort();
444
445 // loop over all muons
446 iterator.Reset();
447 while((candidate = static_cast<Candidate*>(iterator.Next())))
448 {
449 const TLorentzVector &momentum = candidate->Momentum;
450 const TLorentzVector &position = candidate->Position;
451
452
453 pt = momentum.Pt();
454 cosTheta = TMath::Abs(momentum.CosTheta());
455 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
456 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
457 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
458
459 entry = static_cast<Muon*>(branch->NewEntry());
460
461 entry->SetBit(kIsReferenced);
462 entry->SetUniqueID(candidate->GetUniqueID());
463
464 entry->Eta = eta;
465 entry->Phi = momentum.Phi();
466 entry->PT = pt;
467
468 entry->T = position.T()*1.0E-3/c_light;
469
470 // cout<<entry->PT<<","<<entry->T<<endl;
471
472 entry->Charge = candidate->Charge;
473
474 entry->Particle = candidate->GetCandidates()->At(0);
475 }
476}
477
478//------------------------------------------------------------------------------
479
480void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
481{
482 TIter iterator(array);
483 Candidate *candidate = 0, *constituent = 0;
484 Jet *entry = 0;
485 Double_t pt, signPz, cosTheta, eta, rapidity;
486 Double_t ecalEnergy, hcalEnergy;
487 const Double_t c_light = 2.99792458E8;
488
489 array->Sort();
490
491 // loop over all jets
492 iterator.Reset();
493 while((candidate = static_cast<Candidate*>(iterator.Next())))
494 {
495 TIter itConstituents(candidate->GetCandidates());
496
497 const TLorentzVector &momentum = candidate->Momentum;
498 const TLorentzVector &position = candidate->Position;
499
500 pt = momentum.Pt();
501 cosTheta = TMath::Abs(momentum.CosTheta());
502 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
503 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
504 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
505
506 entry = static_cast<Jet*>(branch->NewEntry());
507
508 entry->Eta = eta;
509 entry->Phi = momentum.Phi();
510 entry->PT = pt;
511
512 entry->T = position.T()*1.0E-3/c_light;
513
514 entry->Mass = momentum.M();
515
516 entry->DeltaEta = candidate->DeltaEta;
517 entry->DeltaPhi = candidate->DeltaPhi;
518
519 entry->BTag = candidate->BTag;
520 entry->TauTag = candidate->TauTag;
521
522 entry->Charge = candidate->Charge;
523
524 itConstituents.Reset();
525 entry->Constituents.Clear();
526 ecalEnergy = 0.0;
527 hcalEnergy = 0.0;
528 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
529 {
530 entry->Constituents.Add(constituent);
531 ecalEnergy += constituent->Eem;
532 hcalEnergy += constituent->Ehad;
533 }
534
535 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
536
537 //--- Pile-Up Jet ID variables ----
538
539 entry->NCharged = candidate->NCharged;
540 entry->NNeutrals = candidate->NNeutrals;
541 entry->Beta = candidate->Beta;
542 entry->BetaStar = candidate->BetaStar;
543 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
544 entry->PTD = candidate->PTD;
545 entry->FracPt[0] = candidate->FracPt[0];
546 entry->FracPt[1] = candidate->FracPt[1];
547 entry->FracPt[2] = candidate->FracPt[2];
548 entry->FracPt[3] = candidate->FracPt[3];
549 entry->FracPt[4] = candidate->FracPt[4];
550
551 FillParticles(candidate, &entry->Particles);
552 }
553}
554
555//------------------------------------------------------------------------------
556
557void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
558{
559 Candidate *candidate = 0;
560 MissingET *entry = 0;
561
562 // get the first entry
563 if((candidate = static_cast<Candidate*>(array->At(0))))
564 {
565 const TLorentzVector &momentum = candidate->Momentum;
566
567 entry = static_cast<MissingET*>(branch->NewEntry());
568
569 entry->Eta = (-momentum).Eta();
570 entry->Phi = (-momentum).Phi();
571 entry->MET = momentum.Pt();
572 }
573}
574
575//------------------------------------------------------------------------------
576
577void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
578{
579 Candidate *candidate = 0;
580 ScalarHT *entry = 0;
581
582 // get the first entry
583 if((candidate = static_cast<Candidate*>(array->At(0))))
584 {
585 const TLorentzVector &momentum = candidate->Momentum;
586
587 entry = static_cast<ScalarHT*>(branch->NewEntry());
588
589 entry->HT = momentum.Pt();
590 }
591}
592
593//------------------------------------------------------------------------------
594
595void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
596{
597 TIter iterator(array);
598 Candidate *candidate = 0;
599 Rho *entry = 0;
600
601 // loop over all rho
602 iterator.Reset();
603 while((candidate = static_cast<Candidate*>(iterator.Next())))
604 {
605 const TLorentzVector &momentum = candidate->Momentum;
606
607 entry = static_cast<Rho*>(branch->NewEntry());
608
609 entry->Rho = momentum.E();
610 entry->Edges[0] = candidate->Edges[0];
611 entry->Edges[1] = candidate->Edges[1];
612 }
613}
614
615//------------------------------------------------------------------------------
616
617void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
618{
619 Candidate *candidate = 0;
620 Weight *entry = 0;
621
622 // get the first entry
623 if((candidate = static_cast<Candidate*>(array->At(0))))
624 {
625 const TLorentzVector &momentum = candidate->Momentum;
626
627 entry = static_cast<Weight*>(branch->NewEntry());
628
629 entry->Weight = momentum.E();
630 }
631}
632
633//------------------------------------------------------------------------------
634
635void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
636{
637 TIter iterator(array);
638 Candidate *candidate = 0;
639 HectorHit *entry = 0;
640
641 // loop over all roman pot hits
642 iterator.Reset();
643 while((candidate = static_cast<Candidate*>(iterator.Next())))
644 {
645 const TLorentzVector &position = candidate->Position;
646 const TLorentzVector &momentum = candidate->Momentum;
647
648 entry = static_cast<HectorHit*>(branch->NewEntry());
649
650 entry->E = momentum.E();
651
652 entry->Tx = momentum.Px();
653 entry->Ty = momentum.Py();
654
655 entry->T = position.T();
656
657 entry->X = position.X();
658 entry->Y = position.Y();
659 entry->S = position.Z();
660 }
661}
662
663//------------------------------------------------------------------------------
664
665void TreeWriter::Process()
666{
667 TBranchMap::iterator itBranchMap;
668 ExRootTreeBranch *branch;
669 TProcessMethod method;
670 TObjArray *array;
671
672 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
673 {
674 branch = itBranchMap->first;
675 method = itBranchMap->second.first;
676 array = itBranchMap->second.second;
677
678 (this->*method)(branch, array);
679 }
680}
681
682//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.