Fork me on GitHub

source: git/modules/TreeWriter.cc@ ca55927f

ImprovedOutputFile Timing dual_readout llp 3.0.7
Last change on this file since ca55927f was d7d2da3, checked in by pavel <pavel@…>, 12 years ago

move branches/ModularDelphes to trunk

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