Fork me on GitHub

source: svn/trunk/modules/TreeWriter.cc@ 1120

Last change on this file since 1120 was 1115, checked in by Pavel Demin, 12 years ago

read Rho as TLorentzVector::E()

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