Fork me on GitHub

source: git/modules/TreeWriter.cc@ 2de4dcd

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

fix class map in TreeWriter

  • Property mode set to 100644
File size: 16.0 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
70 TBranchMap::iterator itBranchMap;
71 map< TClass *, TProcessMethod >::iterator itClassMap;
72
73 // read branch configuration and
74 // import array with output from filter/classifier/jetfinder modules
75
76 ExRootConfParam param = GetParam("Branch");
77 Long_t i, size;
78 TString branchName, branchClassName, branchInputArray;
79 TClass *branchClass;
80 TObjArray *array;
81 ExRootTreeBranch *branch;
82
83 size = param.GetSize();
84 for(i = 0; i < size/3; ++i)
85 {
86 branchInputArray = param[i*3].GetString();
87 branchName = param[i*3 + 1].GetString();
88 branchClassName = param[i*3 + 2].GetString();
89
90 branchClass = gROOT->GetClass(branchClassName);
91
92 if(!branchClass)
93 {
94 cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
95 continue;
96 }
97
98 itClassMap = fClassMap.find(branchClass);
99 if(itClassMap == fClassMap.end())
100 {
101 cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
102 continue;
103 }
104
105 array = ImportArray(branchInputArray);
106 branch = NewBranch(branchName, branchClass);
107
108 fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
109 }
110
111}
112
113//------------------------------------------------------------------------------
114
115void TreeWriter::Finish()
116{
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 // loop over all particles
165 iterator.Reset();
166 while((candidate = static_cast<Candidate*>(iterator.Next())))
167 {
168 const TLorentzVector &momentum = candidate->Momentum;
169 const TLorentzVector &position = candidate->Position;
170
171 entry = static_cast<GenParticle*>(branch->NewEntry());
172
173 entry->SetBit(kIsReferenced);
174 entry->SetUniqueID(candidate->GetUniqueID());
175
176 pt = momentum.Pt();
177 cosTheta = TMath::Abs(momentum.CosTheta());
178 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
179 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
180 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
181
182 entry->PID = candidate->PID;
183
184 entry->Status = candidate->Status;
185 entry->IsPU = candidate->IsPU;
186
187 entry->M1 = candidate->M1;
188 entry->M2 = candidate->M2;
189
190 entry->D1 = candidate->D1;
191 entry->D2 = candidate->D2;
192
193 entry->Charge = candidate->Charge;
194 entry->Mass = candidate->Mass;
195
196 entry->E = momentum.E();
197 entry->Px = momentum.Px();
198 entry->Py = momentum.Py();
199 entry->Pz = momentum.Pz();
200
201 entry->Eta = eta;
202 entry->Phi = momentum.Phi();
203 entry->PT = pt;
204
205 entry->Rapidity = rapidity;
206
207 entry->X = position.X();
208 entry->Y = position.Y();
209 entry->Z = position.Z();
210 entry->T = position.T();
211 }
212}
213
214//------------------------------------------------------------------------------
215
216void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
217{
218 TIter iterator(array);
219 Candidate *candidate = 0;
220 Vertex *entry = 0;
221
222 // loop over all vertices
223 iterator.Reset();
224 while((candidate = static_cast<Candidate*>(iterator.Next())))
225 {
226 const TLorentzVector &position = candidate->Position;
227
228 entry = static_cast<Vertex*>(branch->NewEntry());
229
230 entry->X = position.X();
231 entry->Y = position.Y();
232 entry->Z = position.Z();
233 entry->T = position.T();
234 }
235}
236
237//------------------------------------------------------------------------------
238
239void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
240{
241 TIter iterator(array);
242 Candidate *candidate = 0;
243 Candidate *particle = 0;
244 Track *entry = 0;
245 Double_t pt, signz, cosTheta, eta, rapidity;
246
247 // loop over all tracks
248 iterator.Reset();
249 while((candidate = static_cast<Candidate*>(iterator.Next())))
250 {
251 const TLorentzVector &position = candidate->Position;
252
253 cosTheta = TMath::Abs(position.CosTheta());
254 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
255 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
256 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
257
258 entry = static_cast<Track*>(branch->NewEntry());
259
260 entry->SetBit(kIsReferenced);
261 entry->SetUniqueID(candidate->GetUniqueID());
262
263 entry->PID = candidate->PID;
264
265 entry->Charge = candidate->Charge;
266
267 entry->EtaOuter = eta;
268 entry->PhiOuter = position.Phi();
269
270 entry->XOuter = position.X();
271 entry->YOuter = position.Y();
272 entry->ZOuter = position.Z();
273
274 const TLorentzVector &momentum = candidate->Momentum;
275
276 pt = momentum.Pt();
277 cosTheta = TMath::Abs(momentum.CosTheta());
278 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
279 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
280 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
281
282 entry->Eta = eta;
283 entry->Phi = momentum.Phi();
284 entry->PT = pt;
285
286 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
287 const TLorentzVector &initialPosition = particle->Position;
288
289 entry->X = initialPosition.X();
290 entry->Y = initialPosition.Y();
291 entry->Z = initialPosition.Z();
292
293 entry->Particle = particle;
294 }
295}
296
297//------------------------------------------------------------------------------
298
299void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
300{
301 TIter iterator(array);
302 Candidate *candidate = 0;
303 Tower *entry = 0;
304 Double_t pt, signPz, cosTheta, eta, rapidity;
305
306 // loop over all towers
307 iterator.Reset();
308 while((candidate = static_cast<Candidate*>(iterator.Next())))
309 {
310 const TLorentzVector &momentum = candidate->Momentum;
311
312 pt = momentum.Pt();
313 cosTheta = TMath::Abs(momentum.CosTheta());
314 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
315 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
316 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
317
318 entry = static_cast<Tower*>(branch->NewEntry());
319
320 entry->SetBit(kIsReferenced);
321 entry->SetUniqueID(candidate->GetUniqueID());
322
323 entry->Eta = eta;
324 entry->Phi = momentum.Phi();
325 entry->ET = pt;
326 entry->E = momentum.E();
327 entry->Eem = candidate->Eem;
328 entry->Ehad = candidate->Ehad;
329 entry->Edges[0] = candidate->Edges[0];
330 entry->Edges[1] = candidate->Edges[1];
331 entry->Edges[2] = candidate->Edges[2];
332 entry->Edges[3] = candidate->Edges[3];
333
334 FillParticles(candidate, &entry->Particles);
335 }
336}
337
338//------------------------------------------------------------------------------
339
340void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
341{
342 TIter iterator(array);
343 Candidate *candidate = 0;
344 Photon *entry = 0;
345 Double_t pt, signPz, cosTheta, eta, rapidity;
346
347 array->Sort();
348
349 // loop over all photons
350 iterator.Reset();
351 while((candidate = static_cast<Candidate*>(iterator.Next())))
352 {
353 TIter it1(candidate->GetCandidates());
354 const TLorentzVector &momentum = candidate->Momentum;
355
356 pt = momentum.Pt();
357 cosTheta = TMath::Abs(momentum.CosTheta());
358 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
359 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
360 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
361
362 entry = static_cast<Photon*>(branch->NewEntry());
363
364 entry->Eta = eta;
365 entry->Phi = momentum.Phi();
366 entry->PT = pt;
367 entry->E = momentum.E();
368
369 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
370
371 FillParticles(candidate, &entry->Particles);
372 }
373}
374
375//------------------------------------------------------------------------------
376
377void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
378{
379 TIter iterator(array);
380 Candidate *candidate = 0;
381 Electron *entry = 0;
382 Double_t pt, signPz, cosTheta, eta, rapidity;
383
384 array->Sort();
385
386 // loop over all electrons
387 iterator.Reset();
388 while((candidate = static_cast<Candidate*>(iterator.Next())))
389 {
390 const TLorentzVector &momentum = candidate->Momentum;
391
392 pt = momentum.Pt();
393 cosTheta = TMath::Abs(momentum.CosTheta());
394 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
395 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
396 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
397
398 entry = static_cast<Electron*>(branch->NewEntry());
399
400 entry->Eta = eta;
401 entry->Phi = momentum.Phi();
402 entry->PT = pt;
403
404 entry->Charge = candidate->Charge;
405
406 entry->EhadOverEem = 0.0;
407
408 entry->Particle = candidate->GetCandidates()->At(0);
409 }
410}
411
412//------------------------------------------------------------------------------
413
414void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
415{
416 TIter iterator(array);
417 Candidate *candidate = 0;
418 Muon *entry = 0;
419 Double_t pt, signPz, cosTheta, eta, rapidity;
420
421 array->Sort();
422
423 // loop over all muons
424 iterator.Reset();
425 while((candidate = static_cast<Candidate*>(iterator.Next())))
426 {
427 const TLorentzVector &momentum = candidate->Momentum;
428
429 pt = momentum.Pt();
430 cosTheta = TMath::Abs(momentum.CosTheta());
431 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
432 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
433 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
434
435 entry = static_cast<Muon*>(branch->NewEntry());
436
437 entry->SetBit(kIsReferenced);
438 entry->SetUniqueID(candidate->GetUniqueID());
439
440 entry->Eta = eta;
441 entry->Phi = momentum.Phi();
442 entry->PT = pt;
443
444 entry->Charge = candidate->Charge;
445
446 entry->Particle = candidate->GetCandidates()->At(0);
447 }
448}
449
450//------------------------------------------------------------------------------
451
452void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
453{
454 TIter iterator(array);
455 Candidate *candidate = 0, *constituent = 0;
456 Jet *entry = 0;
457 Double_t pt, signPz, cosTheta, eta, rapidity;
458 Double_t ecalEnergy, hcalEnergy;
459
460 array->Sort();
461
462 // loop over all jets
463 iterator.Reset();
464 while((candidate = static_cast<Candidate*>(iterator.Next())))
465 {
466 TIter itConstituents(candidate->GetCandidates());
467 const TLorentzVector &momentum = candidate->Momentum;
468
469 pt = momentum.Pt();
470 cosTheta = TMath::Abs(momentum.CosTheta());
471 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
472 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
473 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
474
475 entry = static_cast<Jet*>(branch->NewEntry());
476
477 entry->Eta = eta;
478 entry->Phi = momentum.Phi();
479 entry->PT = pt;
480
481 entry->Mass = momentum.M();
482
483 entry->DeltaEta = candidate->DeltaEta;
484 entry->DeltaPhi = candidate->DeltaPhi;
485
486 entry->BTag = candidate->BTag;
487 entry->TauTag = candidate->TauTag;
488
489 entry->Charge = candidate->Charge;
490
491 itConstituents.Reset();
492 entry->Constituents.Clear();
493 ecalEnergy = 0.0;
494 hcalEnergy = 0.0;
495 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
496 {
497 entry->Constituents.Add(constituent);
498 ecalEnergy += constituent->Eem;
499 hcalEnergy += constituent->Ehad;
500 }
501
502 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
503
504 FillParticles(candidate, &entry->Particles);
505 }
506}
507
508//------------------------------------------------------------------------------
509
510void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
511{
512 Candidate *candidate = 0;
513 MissingET *entry = 0;
514
515 // get the first entry
516 if((candidate = static_cast<Candidate*>(array->At(0))))
517 {
518 const TLorentzVector &momentum = candidate->Momentum;
519
520 entry = static_cast<MissingET*>(branch->NewEntry());
521
522 entry->Eta = (-momentum).Eta();
523 entry->Phi = (-momentum).Phi();
524 entry->MET = momentum.Pt();
525 }
526}
527
528//------------------------------------------------------------------------------
529
530void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
531{
532 Candidate *candidate = 0;
533 ScalarHT *entry = 0;
534
535 // get the first entry
536 if((candidate = static_cast<Candidate*>(array->At(0))))
537 {
538 const TLorentzVector &momentum = candidate->Momentum;
539
540 entry = static_cast<ScalarHT*>(branch->NewEntry());
541
542 entry->HT = momentum.Pt();
543 }
544}
545
546//------------------------------------------------------------------------------
547
548void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
549{
550 TIter iterator(array);
551 Candidate *candidate = 0;
552 Rho *entry = 0;
553
554 // loop over all rho
555 iterator.Reset();
556 while((candidate = static_cast<Candidate*>(iterator.Next())))
557 {
558 const TLorentzVector &momentum = candidate->Momentum;
559
560 entry = static_cast<Rho*>(branch->NewEntry());
561
562 entry->Rho = momentum.E();
563 entry->Edges[0] = candidate->Edges[0];
564 entry->Edges[1] = candidate->Edges[1];
565 }
566}
567
568//------------------------------------------------------------------------------
569
570void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
571{
572 Candidate *candidate = 0;
573 Weight *entry = 0;
574
575 // get the first entry
576 if((candidate = static_cast<Candidate*>(array->At(0))))
577 {
578 const TLorentzVector &momentum = candidate->Momentum;
579
580 entry = static_cast<Weight*>(branch->NewEntry());
581
582 entry->Weight = momentum.E();
583 }
584}
585
586//------------------------------------------------------------------------------
587
588void TreeWriter::Process()
589{
590 TBranchMap::iterator itBranchMap;
591 ExRootTreeBranch *branch;
592 TProcessMethod method;
593 TObjArray *array;
594
595 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
596 {
597 branch = itBranchMap->first;
598 method = itBranchMap->second.first;
599 array = itBranchMap->second.second;
600
601 (this->*method)(branch, array);
602 }
603}
604
605//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.