Fork me on GitHub

source: git/modules/TreeWriter.cc@ 1150871

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 1150871 was 24d005f, checked in by mselvaggi <mselvaggi@…>, 11 years ago

added PileUpJetID module from Seth Senz

  • Property mode set to 100644
File size: 17.6 KB
RevLine 
[d7d2da3]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;
[2de4dcd]58 fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
[d7d2da3]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;
[71648c2]67 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
[2e229c9]68 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
[d7d2da3]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 }
[22dc7fd]97
[d7d2da3]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
121void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
122{
123 TIter it1(candidate->GetCandidates());
124 it1.Reset();
125 array->Clear();
126 while((candidate = static_cast<Candidate*>(it1.Next())))
127 {
128 TIter it2(candidate->GetCandidates());
129
130 // particle
131 if(candidate->GetCandidates()->GetEntriesFast() == 0)
132 {
133 array->Add(candidate);
134 continue;
135 }
136
137 // track
138 candidate = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
139 if(candidate->GetCandidates()->GetEntriesFast() == 0)
140 {
141 array->Add(candidate);
142 continue;
143 }
144
145 // tower
146 it2.Reset();
147 while((candidate = static_cast<Candidate*>(it2.Next())))
148 {
149 array->Add(candidate->GetCandidates()->At(0));
150 }
151 }
152}
153
154//------------------------------------------------------------------------------
155
156void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
157{
158 TIter iterator(array);
159 Candidate *candidate = 0;
160 GenParticle *entry = 0;
161 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]162
163 const Double_t c_light = 2.99792458E8;
164
[d7d2da3]165 // loop over all particles
166 iterator.Reset();
167 while((candidate = static_cast<Candidate*>(iterator.Next())))
168 {
169 const TLorentzVector &momentum = candidate->Momentum;
170 const TLorentzVector &position = candidate->Position;
171
172 entry = static_cast<GenParticle*>(branch->NewEntry());
173
174 entry->SetBit(kIsReferenced);
175 entry->SetUniqueID(candidate->GetUniqueID());
176
177 pt = momentum.Pt();
178 cosTheta = TMath::Abs(momentum.CosTheta());
179 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
180 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
181 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
182
183 entry->PID = candidate->PID;
184
185 entry->Status = candidate->Status;
186 entry->IsPU = candidate->IsPU;
187
188 entry->M1 = candidate->M1;
189 entry->M2 = candidate->M2;
190
191 entry->D1 = candidate->D1;
192 entry->D2 = candidate->D2;
193
194 entry->Charge = candidate->Charge;
195 entry->Mass = candidate->Mass;
196
197 entry->E = momentum.E();
198 entry->Px = momentum.Px();
199 entry->Py = momentum.Py();
200 entry->Pz = momentum.Pz();
201
202 entry->Eta = eta;
203 entry->Phi = momentum.Phi();
204 entry->PT = pt;
205
206 entry->Rapidity = rapidity;
207
208 entry->X = position.X();
209 entry->Y = position.Y();
210 entry->Z = position.Z();
[22dc7fd]211 entry->T = position.T()*1.0E-3/c_light;
[d7d2da3]212 }
213}
214
215//------------------------------------------------------------------------------
216
[d07e957]217void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
218{
219 TIter iterator(array);
220 Candidate *candidate = 0;
221 Vertex *entry = 0;
[22dc7fd]222
223 const Double_t c_light = 2.99792458E8;
[d07e957]224
225 // loop over all vertices
226 iterator.Reset();
227 while((candidate = static_cast<Candidate*>(iterator.Next())))
228 {
229 const TLorentzVector &position = candidate->Position;
230
231 entry = static_cast<Vertex*>(branch->NewEntry());
232
233 entry->X = position.X();
234 entry->Y = position.Y();
235 entry->Z = position.Z();
[22dc7fd]236 entry->T = position.T()*1.0E-3/c_light;
[d07e957]237 }
238}
239
240//------------------------------------------------------------------------------
241
[d7d2da3]242void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
243{
244 TIter iterator(array);
245 Candidate *candidate = 0;
246 Candidate *particle = 0;
247 Track *entry = 0;
248 Double_t pt, signz, cosTheta, eta, rapidity;
[22dc7fd]249 const Double_t c_light = 2.99792458E8;
250
[d07e957]251 // loop over all tracks
[d7d2da3]252 iterator.Reset();
253 while((candidate = static_cast<Candidate*>(iterator.Next())))
254 {
255 const TLorentzVector &position = candidate->Position;
256
257 cosTheta = TMath::Abs(position.CosTheta());
258 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
259 eta = (cosTheta == 1.0 ? signz*999.9 : position.Eta());
260 rapidity = (cosTheta == 1.0 ? signz*999.9 : position.Rapidity());
261
262 entry = static_cast<Track*>(branch->NewEntry());
263
264 entry->SetBit(kIsReferenced);
265 entry->SetUniqueID(candidate->GetUniqueID());
266
267 entry->PID = candidate->PID;
268
269 entry->Charge = candidate->Charge;
270
271 entry->EtaOuter = eta;
272 entry->PhiOuter = position.Phi();
273
274 entry->XOuter = position.X();
275 entry->YOuter = position.Y();
276 entry->ZOuter = position.Z();
[22dc7fd]277 entry->TOuter = position.T()*1.0E-3/c_light;
[d7d2da3]278
279 const TLorentzVector &momentum = candidate->Momentum;
280
281 pt = momentum.Pt();
282 cosTheta = TMath::Abs(momentum.CosTheta());
283 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
284 eta = (cosTheta == 1.0 ? signz*999.9 : momentum.Eta());
285 rapidity = (cosTheta == 1.0 ? signz*999.9 : momentum.Rapidity());
286
287 entry->Eta = eta;
288 entry->Phi = momentum.Phi();
289 entry->PT = pt;
290
291 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
292 const TLorentzVector &initialPosition = particle->Position;
293
294 entry->X = initialPosition.X();
295 entry->Y = initialPosition.Y();
296 entry->Z = initialPosition.Z();
[22dc7fd]297 entry->T = initialPosition.T()*1.0E-3/c_light;
[d7d2da3]298
299 entry->Particle = particle;
300 }
301}
302
303//------------------------------------------------------------------------------
304
305void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
306{
307 TIter iterator(array);
308 Candidate *candidate = 0;
309 Tower *entry = 0;
310 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]311 const Double_t c_light = 2.99792458E8;
312
[d07e957]313 // loop over all towers
[d7d2da3]314 iterator.Reset();
315 while((candidate = static_cast<Candidate*>(iterator.Next())))
316 {
317 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]318 const TLorentzVector &position = candidate->Position;
319
[d7d2da3]320 pt = momentum.Pt();
321 cosTheta = TMath::Abs(momentum.CosTheta());
322 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
323 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
324 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
325
326 entry = static_cast<Tower*>(branch->NewEntry());
327
328 entry->SetBit(kIsReferenced);
329 entry->SetUniqueID(candidate->GetUniqueID());
330
331 entry->Eta = eta;
332 entry->Phi = momentum.Phi();
333 entry->ET = pt;
334 entry->E = momentum.E();
335 entry->Eem = candidate->Eem;
336 entry->Ehad = candidate->Ehad;
337 entry->Edges[0] = candidate->Edges[0];
338 entry->Edges[1] = candidate->Edges[1];
339 entry->Edges[2] = candidate->Edges[2];
340 entry->Edges[3] = candidate->Edges[3];
[22dc7fd]341
342 entry->T = position.T()*1.0E-3/c_light;
343
[d7d2da3]344 FillParticles(candidate, &entry->Particles);
345 }
346}
347
348//------------------------------------------------------------------------------
349
350void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
351{
352 TIter iterator(array);
353 Candidate *candidate = 0;
354 Photon *entry = 0;
355 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]356 const Double_t c_light = 2.99792458E8;
357
[d7d2da3]358 array->Sort();
359
360 // loop over all photons
361 iterator.Reset();
362 while((candidate = static_cast<Candidate*>(iterator.Next())))
363 {
364 TIter it1(candidate->GetCandidates());
365 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]366 const TLorentzVector &position = candidate->Position;
367
[d7d2da3]368
369 pt = momentum.Pt();
370 cosTheta = TMath::Abs(momentum.CosTheta());
371 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
372 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
373 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
374
375 entry = static_cast<Photon*>(branch->NewEntry());
376
377 entry->Eta = eta;
378 entry->Phi = momentum.Phi();
379 entry->PT = pt;
380 entry->E = momentum.E();
[22dc7fd]381
382 entry->T = position.T()*1.0E-3/c_light;
383
[d7d2da3]384 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9;
385
386 FillParticles(candidate, &entry->Particles);
387 }
388}
389
390//------------------------------------------------------------------------------
391
392void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
393{
394 TIter iterator(array);
395 Candidate *candidate = 0;
396 Electron *entry = 0;
397 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]398 const Double_t c_light = 2.99792458E8;
399
[d7d2da3]400 array->Sort();
401
402 // loop over all electrons
403 iterator.Reset();
404 while((candidate = static_cast<Candidate*>(iterator.Next())))
405 {
406 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]407 const TLorentzVector &position = candidate->Position;
408
[d7d2da3]409 pt = momentum.Pt();
410 cosTheta = TMath::Abs(momentum.CosTheta());
411 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
412 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
413 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
414
415 entry = static_cast<Electron*>(branch->NewEntry());
416
417 entry->Eta = eta;
418 entry->Phi = momentum.Phi();
419 entry->PT = pt;
[22dc7fd]420
421 entry->T = position.T()*1.0E-3/c_light;
422
[d7d2da3]423 entry->Charge = candidate->Charge;
424
425 entry->EhadOverEem = 0.0;
426
427 entry->Particle = candidate->GetCandidates()->At(0);
428 }
429}
430
431//------------------------------------------------------------------------------
432
433void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
434{
435 TIter iterator(array);
436 Candidate *candidate = 0;
437 Muon *entry = 0;
438 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]439
440 const Double_t c_light = 2.99792458E8;
441
[d7d2da3]442 array->Sort();
443
444 // loop over all muons
445 iterator.Reset();
446 while((candidate = static_cast<Candidate*>(iterator.Next())))
447 {
448 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]449 const TLorentzVector &position = candidate->Position;
450
[d7d2da3]451
452 pt = momentum.Pt();
453 cosTheta = TMath::Abs(momentum.CosTheta());
454 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
455 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
456 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
457
458 entry = static_cast<Muon*>(branch->NewEntry());
459
460 entry->SetBit(kIsReferenced);
461 entry->SetUniqueID(candidate->GetUniqueID());
462
463 entry->Eta = eta;
464 entry->Phi = momentum.Phi();
465 entry->PT = pt;
466
[22dc7fd]467 entry->T = position.T()*1.0E-3/c_light;
468
469 // cout<<entry->PT<<","<<entry->T<<endl;
470
[d7d2da3]471 entry->Charge = candidate->Charge;
472
473 entry->Particle = candidate->GetCandidates()->At(0);
474 }
475}
476
477//------------------------------------------------------------------------------
478
479void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
480{
481 TIter iterator(array);
482 Candidate *candidate = 0, *constituent = 0;
483 Jet *entry = 0;
484 Double_t pt, signPz, cosTheta, eta, rapidity;
485 Double_t ecalEnergy, hcalEnergy;
[22dc7fd]486 const Double_t c_light = 2.99792458E8;
487
[d7d2da3]488 array->Sort();
489
490 // loop over all jets
491 iterator.Reset();
492 while((candidate = static_cast<Candidate*>(iterator.Next())))
493 {
494 TIter itConstituents(candidate->GetCandidates());
[22dc7fd]495
[d7d2da3]496 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]497 const TLorentzVector &position = candidate->Position;
498
[d7d2da3]499 pt = momentum.Pt();
500 cosTheta = TMath::Abs(momentum.CosTheta());
501 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
502 eta = (cosTheta == 1.0 ? signPz*999.9 : momentum.Eta());
503 rapidity = (cosTheta == 1.0 ? signPz*999.9 : momentum.Rapidity());
504
505 entry = static_cast<Jet*>(branch->NewEntry());
506
507 entry->Eta = eta;
508 entry->Phi = momentum.Phi();
509 entry->PT = pt;
510
[22dc7fd]511 entry->T = position.T()*1.0E-3/c_light;
512
[d7d2da3]513 entry->Mass = momentum.M();
514
515 entry->DeltaEta = candidate->DeltaEta;
516 entry->DeltaPhi = candidate->DeltaPhi;
517
518 entry->BTag = candidate->BTag;
519 entry->TauTag = candidate->TauTag;
520
521 entry->Charge = candidate->Charge;
522
523 itConstituents.Reset();
524 entry->Constituents.Clear();
525 ecalEnergy = 0.0;
526 hcalEnergy = 0.0;
527 while((constituent = static_cast<Candidate*>(itConstituents.Next())))
528 {
529 entry->Constituents.Add(constituent);
530 ecalEnergy += constituent->Eem;
531 hcalEnergy += constituent->Ehad;
532 }
533
534 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9;
[24d005f]535
536 //--- Pile-Up Jet ID variables ----
537
538 entry->NCharged = candidate->NCharged;
539 entry->NNeutrals = candidate->NNeutrals;
540 entry->Beta = candidate->Beta;
541 entry->BetaStar = candidate->BetaStar;
542 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
543 entry->PTD = candidate->PTD;
544 entry->FracPt[0] = candidate->FracPt[0];
545 entry->FracPt[1] = candidate->FracPt[1];
546 entry->FracPt[2] = candidate->FracPt[2];
547 entry->FracPt[3] = candidate->FracPt[3];
548 entry->FracPt[4] = candidate->FracPt[4];
[d7d2da3]549
550 FillParticles(candidate, &entry->Particles);
551 }
552}
553
554//------------------------------------------------------------------------------
555
556void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
557{
558 Candidate *candidate = 0;
559 MissingET *entry = 0;
560
561 // get the first entry
562 if((candidate = static_cast<Candidate*>(array->At(0))))
563 {
564 const TLorentzVector &momentum = candidate->Momentum;
565
566 entry = static_cast<MissingET*>(branch->NewEntry());
567
[3b465ca]568 entry->Eta = (-momentum).Eta();
[d7d2da3]569 entry->Phi = (-momentum).Phi();
570 entry->MET = momentum.Pt();
571 }
572}
573
574//------------------------------------------------------------------------------
575
576void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
577{
578 Candidate *candidate = 0;
579 ScalarHT *entry = 0;
580
581 // get the first entry
582 if((candidate = static_cast<Candidate*>(array->At(0))))
583 {
584 const TLorentzVector &momentum = candidate->Momentum;
585
586 entry = static_cast<ScalarHT*>(branch->NewEntry());
587
588 entry->HT = momentum.Pt();
589 }
590}
591
592//------------------------------------------------------------------------------
593
[71648c2]594void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
595{
[3b465ca]596 TIter iterator(array);
[71648c2]597 Candidate *candidate = 0;
598 Rho *entry = 0;
599
[3b465ca]600 // loop over all rho
601 iterator.Reset();
602 while((candidate = static_cast<Candidate*>(iterator.Next())))
[71648c2]603 {
604 const TLorentzVector &momentum = candidate->Momentum;
605
606 entry = static_cast<Rho*>(branch->NewEntry());
607
[8608ef8]608 entry->Rho = momentum.E();
[3b465ca]609 entry->Edges[0] = candidate->Edges[0];
610 entry->Edges[1] = candidate->Edges[1];
[71648c2]611 }
612}
613
614//------------------------------------------------------------------------------
615
[2e229c9]616void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
617{
618 Candidate *candidate = 0;
619 Weight *entry = 0;
620
621 // get the first entry
622 if((candidate = static_cast<Candidate*>(array->At(0))))
623 {
624 const TLorentzVector &momentum = candidate->Momentum;
625
626 entry = static_cast<Weight*>(branch->NewEntry());
627
628 entry->Weight = momentum.E();
629 }
630}
631
632//------------------------------------------------------------------------------
633
[d7d2da3]634void TreeWriter::Process()
635{
636 TBranchMap::iterator itBranchMap;
637 ExRootTreeBranch *branch;
638 TProcessMethod method;
639 TObjArray *array;
640
641 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
642 {
643 branch = itBranchMap->first;
644 method = itBranchMap->second.first;
645 array = itBranchMap->second.second;
646
647 (this->*method)(branch, array);
648 }
649}
650
651//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.