Fork me on GitHub

source: git/modules/TreeWriter.cc@ 78e55cf

Last change on this file since 78e55cf was 0060caa, checked in by michele <michele.selvaggi@…>, 4 years ago

initialPos in PF and Tracks given by cand, not particle pointer

  • Property mode set to 100644
File size: 28.9 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
[1fa50c2]9 *
[b443089]10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
[1fa50c2]14 *
[b443089]15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
[d7d2da3]19/** \class TreeWriter
20 *
21 * Fills ROOT tree branches.
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/TreeWriter.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
[341014c]34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
[d7d2da3]36#include "ExRootAnalysis/ExRootTreeBranch.h"
37
38#include "TDatabasePDG.h"
[341014c]39#include "TFormula.h"
[d7d2da3]40#include "TLorentzVector.h"
[341014c]41#include "TMath.h"
42#include "TObjArray.h"
43#include "TROOT.h"
44#include "TRandom3.h"
45#include "TString.h"
[d7d2da3]46
47#include <algorithm>
48#include <iostream>
49#include <sstream>
[341014c]50#include <stdexcept>
[d7d2da3]51
52using namespace std;
53
54//------------------------------------------------------------------------------
55
56TreeWriter::TreeWriter()
57{
58}
59
60//------------------------------------------------------------------------------
61
62TreeWriter::~TreeWriter()
63{
64}
65
66//------------------------------------------------------------------------------
67
68void TreeWriter::Init()
69{
70 fClassMap[GenParticle::Class()] = &TreeWriter::ProcessParticles;
[2de4dcd]71 fClassMap[Vertex::Class()] = &TreeWriter::ProcessVertices;
[d7d2da3]72 fClassMap[Track::Class()] = &TreeWriter::ProcessTracks;
73 fClassMap[Tower::Class()] = &TreeWriter::ProcessTowers;
[4d7014e]74 fClassMap[ParticleFlowCandidate::Class()] = &TreeWriter::ProcessParticleFlowCandidates;
[d7d2da3]75 fClassMap[Photon::Class()] = &TreeWriter::ProcessPhotons;
76 fClassMap[Electron::Class()] = &TreeWriter::ProcessElectrons;
77 fClassMap[Muon::Class()] = &TreeWriter::ProcessMuons;
78 fClassMap[Jet::Class()] = &TreeWriter::ProcessJets;
79 fClassMap[MissingET::Class()] = &TreeWriter::ProcessMissingET;
80 fClassMap[ScalarHT::Class()] = &TreeWriter::ProcessScalarHT;
[71648c2]81 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho;
[2e229c9]82 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight;
[8f7db23]83 fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit;
[d7d2da3]84
85 TBranchMap::iterator itBranchMap;
[341014c]86 map<TClass *, TProcessMethod>::iterator itClassMap;
[d7d2da3]87
88 // read branch configuration and
89 // import array with output from filter/classifier/jetfinder modules
90
91 ExRootConfParam param = GetParam("Branch");
92 Long_t i, size;
93 TString branchName, branchClassName, branchInputArray;
94 TClass *branchClass;
95 TObjArray *array;
96 ExRootTreeBranch *branch;
97
98 size = param.GetSize();
[341014c]99 for(i = 0; i < size / 3; ++i)
[d7d2da3]100 {
[341014c]101 branchInputArray = param[i * 3].GetString();
102 branchName = param[i * 3 + 1].GetString();
103 branchClassName = param[i * 3 + 2].GetString();
[d7d2da3]104
105 branchClass = gROOT->GetClass(branchClassName);
106
107 if(!branchClass)
108 {
109 cout << "** ERROR: cannot find class '" << branchClassName << "'" << endl;
110 continue;
111 }
[8f7db23]112
[d7d2da3]113 itClassMap = fClassMap.find(branchClass);
114 if(itClassMap == fClassMap.end())
115 {
116 cout << "** ERROR: cannot create branch for class '" << branchClassName << "'" << endl;
117 continue;
118 }
119
120 array = ImportArray(branchInputArray);
121 branch = NewBranch(branchName, branchClass);
122
123 fBranchMap.insert(make_pair(branch, make_pair(itClassMap->second, array)));
124 }
[dc883b4]125
126 param = GetParam("Info");
127 TString infoName;
128 Double_t infoValue;
129
130 size = param.GetSize();
131 for(i = 0; i < size / 2; ++i)
132 {
133 infoName = param[i * 2].GetString();
134 infoValue = param[i * 2 + 1].GetDouble();
135
136 AddInfo(infoName, infoValue);
137 }
[d7d2da3]138}
139
140//------------------------------------------------------------------------------
141
142void TreeWriter::Finish()
143{
144}
145
146//------------------------------------------------------------------------------
147
148void TreeWriter::FillParticles(Candidate *candidate, TRefArray *array)
149{
150 TIter it1(candidate->GetCandidates());
151 it1.Reset();
152 array->Clear();
[dc883b4]153
[341014c]154 while((candidate = static_cast<Candidate *>(it1.Next())))
[d7d2da3]155 {
156 TIter it2(candidate->GetCandidates());
157
158 // particle
159 if(candidate->GetCandidates()->GetEntriesFast() == 0)
160 {
161 array->Add(candidate);
162 continue;
163 }
164
165 // track
[341014c]166 candidate = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
[d7d2da3]167 if(candidate->GetCandidates()->GetEntriesFast() == 0)
168 {
169 array->Add(candidate);
170 continue;
171 }
172
173 // tower
174 it2.Reset();
[341014c]175 while((candidate = static_cast<Candidate *>(it2.Next())))
[d7d2da3]176 {
177 array->Add(candidate->GetCandidates()->At(0));
178 }
179 }
180}
181
182//------------------------------------------------------------------------------
183
184void TreeWriter::ProcessParticles(ExRootTreeBranch *branch, TObjArray *array)
185{
186 TIter iterator(array);
187 Candidate *candidate = 0;
188 GenParticle *entry = 0;
189 Double_t pt, signPz, cosTheta, eta, rapidity;
[8f7db23]190
[22dc7fd]191 const Double_t c_light = 2.99792458E8;
[8f7db23]192
[d7d2da3]193 // loop over all particles
194 iterator.Reset();
[341014c]195 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]196 {
197 const TLorentzVector &momentum = candidate->Momentum;
198 const TLorentzVector &position = candidate->Position;
199
[341014c]200 entry = static_cast<GenParticle *>(branch->NewEntry());
[d7d2da3]201
202 entry->SetBit(kIsReferenced);
203 entry->SetUniqueID(candidate->GetUniqueID());
204
205 pt = momentum.Pt();
206 cosTheta = TMath::Abs(momentum.CosTheta());
207 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]208 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
209 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]210
211 entry->PID = candidate->PID;
212
213 entry->Status = candidate->Status;
214 entry->IsPU = candidate->IsPU;
215
216 entry->M1 = candidate->M1;
217 entry->M2 = candidate->M2;
218
219 entry->D1 = candidate->D1;
220 entry->D2 = candidate->D2;
221
222 entry->Charge = candidate->Charge;
223 entry->Mass = candidate->Mass;
224
225 entry->E = momentum.E();
226 entry->Px = momentum.Px();
227 entry->Py = momentum.Py();
228 entry->Pz = momentum.Pz();
229
230 entry->Eta = eta;
231 entry->Phi = momentum.Phi();
232 entry->PT = pt;
233
234 entry->Rapidity = rapidity;
235
236 entry->X = position.X();
237 entry->Y = position.Y();
238 entry->Z = position.Z();
[341014c]239 entry->T = position.T() * 1.0E-3 / c_light;
[d7d2da3]240 }
241}
242
243//------------------------------------------------------------------------------
244
[d07e957]245void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array)
246{
247 TIter iterator(array);
[5496767]248 Candidate *candidate = 0, *constituent = 0;
[d07e957]249 Vertex *entry = 0;
[8f7db23]250
[22dc7fd]251 const Double_t c_light = 2.99792458E8;
[5496767]252
[332025f]253 Double_t x, y, z, t, xError, yError, zError, tError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2;
[2600216]254 UInt_t index, ndf;
[5496767]255
[3e2bb2b]256 CompBase *compare = Candidate::fgCompare;
257 Candidate::fgCompare = CompSumPT2<Candidate>::Instance();
[3c46e17]258 array->Sort();
[3e2bb2b]259 Candidate::fgCompare = compare;
[5496767]260
[d07e957]261 // loop over all vertices
262 iterator.Reset();
[341014c]263 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d07e957]264 {
[5496767]265
[2600216]266 index = candidate->ClusterIndex;
267 ndf = candidate->ClusterNDF;
268 sigma = candidate->ClusterSigma;
269 sumPT2 = candidate->SumPT2;
270 btvSumPT2 = candidate->BTVSumPT2;
271 genDeltaZ = candidate->GenDeltaZ;
272 genSumPT2 = candidate->GenSumPT2;
[5496767]273
[7bcca65]274 x = candidate->Position.X();
275 y = candidate->Position.Y();
276 z = candidate->Position.Z();
[341014c]277 t = candidate->Position.T() * 1.0E-3 / c_light;
[5496767]278
[341014c]279 xError = candidate->PositionError.X();
280 yError = candidate->PositionError.Y();
281 zError = candidate->PositionError.Z();
282 tError = candidate->PositionError.T() * 1.0E-3 / c_light;
[d07e957]283
[341014c]284 entry = static_cast<Vertex *>(branch->NewEntry());
[d07e957]285
[2600216]286 entry->Index = index;
287 entry->NDF = ndf;
288 entry->Sigma = sigma;
289 entry->SumPT2 = sumPT2;
290 entry->BTVSumPT2 = btvSumPT2;
291 entry->GenDeltaZ = genDeltaZ;
292 entry->GenSumPT2 = genSumPT2;
[5496767]293
[2600216]294 entry->X = x;
295 entry->Y = y;
296 entry->Z = z;
297 entry->T = t;
[5496767]298
[2600216]299 entry->ErrorX = xError;
300 entry->ErrorY = yError;
301 entry->ErrorZ = zError;
[332025f]302 entry->ErrorT = tError;
[5496767]303
[3e2bb2b]304 TIter itConstituents(candidate->GetCandidates());
305 itConstituents.Reset();
306 entry->Constituents.Clear();
[341014c]307 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
[5496767]308 {
309 entry->Constituents.Add(constituent);
310 }
[d07e957]311 }
312}
313
314//------------------------------------------------------------------------------
315
[d7d2da3]316void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array)
317{
318 TIter iterator(array);
319 Candidate *candidate = 0;
320 Candidate *particle = 0;
321 Track *entry = 0;
[fd4b326]322 Double_t pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
[22dc7fd]323 const Double_t c_light = 2.99792458E8;
[8f7db23]324
[d07e957]325 // loop over all tracks
[d7d2da3]326 iterator.Reset();
[341014c]327 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]328 {
329 const TLorentzVector &position = candidate->Position;
330
331 cosTheta = TMath::Abs(position.CosTheta());
332 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]333 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
334 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
[d7d2da3]335
[341014c]336 entry = static_cast<Track *>(branch->NewEntry());
[d7d2da3]337
338 entry->SetBit(kIsReferenced);
339 entry->SetUniqueID(candidate->GetUniqueID());
340
341 entry->PID = candidate->PID;
342
343 entry->Charge = candidate->Charge;
344
345 entry->EtaOuter = eta;
346 entry->PhiOuter = position.Phi();
347
348 entry->XOuter = position.X();
349 entry->YOuter = position.Y();
350 entry->ZOuter = position.Z();
[341014c]351 entry->TOuter = position.T() * 1.0E-3 / c_light;
[5496767]352
[acd0621]353 entry->L = candidate->L;
[5496767]354
[341014c]355 entry->D0 = candidate->D0;
356 entry->DZ = candidate->DZ;
[a95da74]357 entry->Nclusters = candidate->Nclusters;
[2a3eb22]358
[341014c]359 entry->ErrorP = candidate->ErrorP;
360 entry->ErrorPT = candidate->ErrorPT;
[2671df6]361
362 // diagonal covariance matrix terms
363 entry->ErrorD0 = candidate->ErrorD0;
364 entry->ErrorC = candidate->ErrorC;
[341014c]365 entry->ErrorPhi = candidate->ErrorPhi;
[2671df6]366 entry->ErrorDZ = candidate->ErrorDZ;
367 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
368
369 // add some offdiagonal covariance matrix elements
370 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
371 entry->ErrorD0C = candidate->TrackCovariance(0,2);
372 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
373 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
374 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
375 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
376 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
377 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
378 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
379 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
[5496767]380
[e4c3fef]381 entry->Xd = candidate->Xd;
382 entry->Yd = candidate->Yd;
383 entry->Zd = candidate->Zd;
[9040259]384
[d7d2da3]385 const TLorentzVector &momentum = candidate->Momentum;
386
387 pt = momentum.Pt();
[2a3eb22]388 p = momentum.P();
389 phi = momentum.Phi();
[fd4b326]390 m = momentum.M();
[341014c]391 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
[2a3eb22]392
[d7d2da3]393 cosTheta = TMath::Abs(momentum.CosTheta());
394 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]395 eta = (cosTheta == 1.0 ? signz * 999.9 : momentum.Eta());
396 rapidity = (cosTheta == 1.0 ? signz * 999.9 : momentum.Rapidity());
[d7d2da3]397
[c1ea422]398 entry->P = p;
[341014c]399 entry->PT = pt;
[d7d2da3]400 entry->Eta = eta;
[2a3eb22]401 entry->Phi = phi;
402 entry->CtgTheta = ctgTheta;
[17cd992]403 entry->C = candidate->C;
[fd4b326]404 entry->Mass = m;
[5496767]405
[0060caa]406 //particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
407 //const TLorentzVector &initialPosition = particle->Position;
408 const TLorentzVector &initialPosition = candidate->InitialPosition;
[d7d2da3]409
410 entry->X = initialPosition.X();
411 entry->Y = initialPosition.Y();
412 entry->Z = initialPosition.Z();
[341014c]413 entry->T = initialPosition.T() * 1.0E-3 / c_light;
[d7d2da3]414
415 entry->Particle = particle;
[2600216]416
417 entry->VertexIndex = candidate->ClusterIndex;
[d7d2da3]418 }
419}
420
421//------------------------------------------------------------------------------
422
423void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
424{
425 TIter iterator(array);
426 Candidate *candidate = 0;
427 Tower *entry = 0;
428 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]429 const Double_t c_light = 2.99792458E8;
[8f7db23]430
[d07e957]431 // loop over all towers
[d7d2da3]432 iterator.Reset();
[341014c]433 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]434 {
435 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]436 const TLorentzVector &position = candidate->Position;
[8f7db23]437
[d7d2da3]438 pt = momentum.Pt();
439 cosTheta = TMath::Abs(momentum.CosTheta());
440 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]441 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
442 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]443
[341014c]444 entry = static_cast<Tower *>(branch->NewEntry());
[d7d2da3]445
446 entry->SetBit(kIsReferenced);
447 entry->SetUniqueID(candidate->GetUniqueID());
448
449 entry->Eta = eta;
450 entry->Phi = momentum.Phi();
451 entry->ET = pt;
452 entry->E = momentum.E();
453 entry->Eem = candidate->Eem;
454 entry->Ehad = candidate->Ehad;
455 entry->Edges[0] = candidate->Edges[0];
456 entry->Edges[1] = candidate->Edges[1];
457 entry->Edges[2] = candidate->Edges[2];
458 entry->Edges[3] = candidate->Edges[3];
[8f7db23]459
[341014c]460 entry->T = position.T() * 1.0E-3 / c_light;
[839deb7]461 entry->NTimeHits = candidate->NTimeHits;
[8f7db23]462
[d7d2da3]463 FillParticles(candidate, &entry->Particles);
464 }
465}
466
467//------------------------------------------------------------------------------
468
[4d7014e]469void TreeWriter::ProcessParticleFlowCandidates(ExRootTreeBranch *branch, TObjArray *array)
470{
471
472 TIter iterator(array);
473 Candidate *candidate = 0;
474 Candidate *particle = 0;
475 ParticleFlowCandidate *entry = 0;
[fd4b326]476 Double_t e, pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
[4d7014e]477 const Double_t c_light = 2.99792458E8;
478
479 // loop over all tracks
480 iterator.Reset();
481 while((candidate = static_cast<Candidate *>(iterator.Next())))
482 {
483 const TLorentzVector &position = candidate->Position;
484
485 cosTheta = TMath::Abs(position.CosTheta());
486 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
487 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
488 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
489
490 entry = static_cast<ParticleFlowCandidate *>(branch->NewEntry());
491
492 entry->SetBit(kIsReferenced);
493 entry->SetUniqueID(candidate->GetUniqueID());
494
495 entry->PID = candidate->PID;
496
497 entry->Charge = candidate->Charge;
498
499 entry->EtaOuter = eta;
500 entry->PhiOuter = position.Phi();
501
502 entry->XOuter = position.X();
503 entry->YOuter = position.Y();
504 entry->ZOuter = position.Z();
505 entry->TOuter = position.T() * 1.0E-3 / c_light;
506
507 entry->L = candidate->L;
508
509 entry->D0 = candidate->D0;
510 entry->DZ = candidate->DZ;
[a95da74]511 entry->Nclusters = candidate->Nclusters;
[4d7014e]512
513 entry->ErrorP = candidate->ErrorP;
514 entry->ErrorPT = candidate->ErrorPT;
515 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
[2671df6]516
517
518 // diagonal covariance matrix terms
519
520 entry->ErrorD0 = candidate->ErrorD0;
521 entry->ErrorC = candidate->ErrorC;
[4d7014e]522 entry->ErrorPhi = candidate->ErrorPhi;
[2671df6]523 entry->ErrorDZ = candidate->ErrorDZ;
524 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
525
526 // add some offdiagonal covariance matrix elements
527 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
528 entry->ErrorD0C = candidate->TrackCovariance(0,2);
529 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
530 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
531 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
532 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
533 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
534 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
535 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
536 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
537
[4d7014e]538 entry->Xd = candidate->Xd;
539 entry->Yd = candidate->Yd;
540 entry->Zd = candidate->Zd;
541
542 const TLorentzVector &momentum = candidate->Momentum;
543
544 e = momentum.E();
545 pt = momentum.Pt();
546 p = momentum.P();
547 phi = momentum.Phi();
[fd4b326]548 m = momentum.M();
[4d7014e]549 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
550
551 entry->E = e;
552 entry->P = p;
553 entry->PT = pt;
554 entry->Eta = eta;
555 entry->Phi = phi;
556 entry->CtgTheta = ctgTheta;
[17cd992]557 entry->C = candidate->C;
[fd4b326]558 entry->Mass = m;
[dc883b4]559
[0060caa]560 //particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
561 //const TLorentzVector &initialPosition = particle->Position;
562 const TLorentzVector &initialPosition = candidate->InitialPosition;
[4d7014e]563
564 entry->X = initialPosition.X();
565 entry->Y = initialPosition.Y();
566 entry->Z = initialPosition.Z();
567 entry->T = initialPosition.T() * 1.0E-3 / c_light;
568
569 entry->VertexIndex = candidate->ClusterIndex;
570
571 entry->Eem = candidate->Eem;
572 entry->Ehad = candidate->Ehad;
573 entry->Edges[0] = candidate->Edges[0];
574 entry->Edges[1] = candidate->Edges[1];
575 entry->Edges[2] = candidate->Edges[2];
576 entry->Edges[3] = candidate->Edges[3];
577
[0060caa]578 //entry->T = position.T() * 1.0E-3 / c_light;
[4d7014e]579 entry->NTimeHits = candidate->NTimeHits;
580
581 FillParticles(candidate, &entry->Particles);
582
583 }
584}
585
586//------------------------------------------------------------------------------
587
[d7d2da3]588void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
589{
590 TIter iterator(array);
591 Candidate *candidate = 0;
592 Photon *entry = 0;
593 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]594 const Double_t c_light = 2.99792458E8;
[8f7db23]595
[d7d2da3]596 array->Sort();
597
598 // loop over all photons
599 iterator.Reset();
[341014c]600 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]601 {
602 TIter it1(candidate->GetCandidates());
603 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]604 const TLorentzVector &position = candidate->Position;
[8f7db23]605
[d7d2da3]606 pt = momentum.Pt();
607 cosTheta = TMath::Abs(momentum.CosTheta());
608 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]609 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
610 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]611
[341014c]612 entry = static_cast<Photon *>(branch->NewEntry());
[d7d2da3]613
614 entry->Eta = eta;
615 entry->Phi = momentum.Phi();
616 entry->PT = pt;
617 entry->E = momentum.E();
[341014c]618 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]619
[d5cae2b]620 // Isolation variables
[9040259]621
[d5cae2b]622 entry->IsolationVar = candidate->IsolationVar;
[341014c]623 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
624 entry->SumPtCharged = candidate->SumPtCharged;
625 entry->SumPtNeutral = candidate->SumPtNeutral;
626 entry->SumPtChargedPU = candidate->SumPtChargedPU;
627 entry->SumPt = candidate->SumPt;
[d5cae2b]628
[341014c]629 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad / candidate->Eem : 999.9;
[d7d2da3]630
[0e0f211]631 // 1: prompt -- 2: non prompt -- 3: fake
632 entry->Status = candidate->Status;
633
[d7d2da3]634 FillParticles(candidate, &entry->Particles);
635 }
636}
637
638//------------------------------------------------------------------------------
639
640void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
641{
642 TIter iterator(array);
643 Candidate *candidate = 0;
644 Electron *entry = 0;
645 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]646 const Double_t c_light = 2.99792458E8;
[8f7db23]647
[d7d2da3]648 array->Sort();
649
650 // loop over all electrons
651 iterator.Reset();
[341014c]652 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]653 {
654 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]655 const TLorentzVector &position = candidate->Position;
[8f7db23]656
[d7d2da3]657 pt = momentum.Pt();
658 cosTheta = TMath::Abs(momentum.CosTheta());
659 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]660 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
661 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]662
[341014c]663 entry = static_cast<Electron *>(branch->NewEntry());
[d7d2da3]664
665 entry->Eta = eta;
666 entry->Phi = momentum.Phi();
667 entry->PT = pt;
[8f7db23]668
[341014c]669 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]670
[0518688]671 // displacement
[341014c]672 entry->D0 = candidate->D0;
673 entry->ErrorD0 = candidate->ErrorD0;
674 entry->DZ = candidate->DZ;
675 entry->ErrorDZ = candidate->ErrorDZ;
[9040259]676
[0518688]677 // Isolation variables
[d5cae2b]678 entry->IsolationVar = candidate->IsolationVar;
[341014c]679 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
680 entry->SumPtCharged = candidate->SumPtCharged;
681 entry->SumPtNeutral = candidate->SumPtNeutral;
682 entry->SumPtChargedPU = candidate->SumPtChargedPU;
683 entry->SumPt = candidate->SumPt;
[d5cae2b]684
[d7d2da3]685 entry->Charge = candidate->Charge;
686
687 entry->EhadOverEem = 0.0;
688
689 entry->Particle = candidate->GetCandidates()->At(0);
690 }
691}
692
693//------------------------------------------------------------------------------
694
695void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
696{
697 TIter iterator(array);
698 Candidate *candidate = 0;
699 Muon *entry = 0;
700 Double_t pt, signPz, cosTheta, eta, rapidity;
[8f7db23]701
[22dc7fd]702 const Double_t c_light = 2.99792458E8;
[8f7db23]703
[d7d2da3]704 array->Sort();
705
706 // loop over all muons
707 iterator.Reset();
[341014c]708 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]709 {
710 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]711 const TLorentzVector &position = candidate->Position;
[8f7db23]712
[d7d2da3]713 pt = momentum.Pt();
714 cosTheta = TMath::Abs(momentum.CosTheta());
715 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]716 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
717 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]718
[341014c]719 entry = static_cast<Muon *>(branch->NewEntry());
[d7d2da3]720
721 entry->SetBit(kIsReferenced);
722 entry->SetUniqueID(candidate->GetUniqueID());
723
724 entry->Eta = eta;
725 entry->Phi = momentum.Phi();
726 entry->PT = pt;
727
[341014c]728 entry->T = position.T() * 1.0E-3 / c_light;
[0518688]729
730 // displacement
[341014c]731 entry->D0 = candidate->D0;
732 entry->ErrorD0 = candidate->ErrorD0;
733 entry->DZ = candidate->DZ;
734 entry->ErrorDZ = candidate->ErrorDZ;
[0518688]735
[d5cae2b]736 // Isolation variables
[9040259]737
[d5cae2b]738 entry->IsolationVar = candidate->IsolationVar;
[341014c]739 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
740 entry->SumPtCharged = candidate->SumPtCharged;
741 entry->SumPtNeutral = candidate->SumPtNeutral;
742 entry->SumPtChargedPU = candidate->SumPtChargedPU;
743 entry->SumPt = candidate->SumPt;
[8f7db23]744
[d7d2da3]745 entry->Charge = candidate->Charge;
746
747 entry->Particle = candidate->GetCandidates()->At(0);
748 }
749}
750
751//------------------------------------------------------------------------------
752
753void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
754{
755 TIter iterator(array);
756 Candidate *candidate = 0, *constituent = 0;
757 Jet *entry = 0;
758 Double_t pt, signPz, cosTheta, eta, rapidity;
759 Double_t ecalEnergy, hcalEnergy;
[22dc7fd]760 const Double_t c_light = 2.99792458E8;
[de6d698]761 Int_t i;
[8f7db23]762
[d7d2da3]763 array->Sort();
764
765 // loop over all jets
766 iterator.Reset();
[341014c]767 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]768 {
769 TIter itConstituents(candidate->GetCandidates());
[8f7db23]770
[d7d2da3]771 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]772 const TLorentzVector &position = candidate->Position;
[8f7db23]773
[d7d2da3]774 pt = momentum.Pt();
775 cosTheta = TMath::Abs(momentum.CosTheta());
776 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]777 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
778 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]779
[341014c]780 entry = static_cast<Jet *>(branch->NewEntry());
[d7d2da3]781
782 entry->Eta = eta;
783 entry->Phi = momentum.Phi();
784 entry->PT = pt;
785
[341014c]786 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]787
[d7d2da3]788 entry->Mass = momentum.M();
789
[ba1f1ee]790 entry->Area = candidate->Area;
791
[d7d2da3]792 entry->DeltaEta = candidate->DeltaEta;
793 entry->DeltaPhi = candidate->DeltaPhi;
794
[fe0273c]795 entry->Flavor = candidate->Flavor;
796 entry->FlavorAlgo = candidate->FlavorAlgo;
797 entry->FlavorPhys = candidate->FlavorPhys;
798
[d7d2da3]799 entry->BTag = candidate->BTag;
[9040259]800
801 entry->BTagAlgo = candidate->BTagAlgo;
[fe0273c]802 entry->BTagPhys = candidate->BTagPhys;
[9040259]803
[d7d2da3]804 entry->TauTag = candidate->TauTag;
[7429c6a]805 entry->TauWeight = candidate->TauWeight;
[d7d2da3]806
807 entry->Charge = candidate->Charge;
808
809 itConstituents.Reset();
810 entry->Constituents.Clear();
811 ecalEnergy = 0.0;
812 hcalEnergy = 0.0;
[341014c]813 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
[d7d2da3]814 {
815 entry->Constituents.Add(constituent);
816 ecalEnergy += constituent->Eem;
817 hcalEnergy += constituent->Ehad;
818 }
819
[341014c]820 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy / ecalEnergy : 999.9;
[8f7db23]821
[24d005f]822 //--- Pile-Up Jet ID variables ----
823
824 entry->NCharged = candidate->NCharged;
825 entry->NNeutrals = candidate->NNeutrals;
[c614dd7]826
827 entry->NeutralEnergyFraction = candidate->NeutralEnergyFraction;
828 entry->ChargedEnergyFraction = candidate->ChargedEnergyFraction;
[24d005f]829 entry->Beta = candidate->Beta;
830 entry->BetaStar = candidate->BetaStar;
831 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
832 entry->PTD = candidate->PTD;
[9040259]833
[de6d698]834 //--- Sub-structure variables ----
[9040259]835
836 entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
837 entry->NSubJetsPruned = candidate->NSubJetsPruned;
[de6d698]838 entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
[9040259]839
[341014c]840 entry->SoftDroppedJet = candidate->SoftDroppedJet;
841 entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1;
[ba75867]842 entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
843
[de6d698]844 for(i = 0; i < 5; i++)
845 {
[341014c]846 entry->FracPt[i] = candidate->FracPt[i];
847 entry->Tau[i] = candidate->Tau[i];
848 entry->TrimmedP4[i] = candidate->TrimmedP4[i];
849 entry->PrunedP4[i] = candidate->PrunedP4[i];
850 entry->SoftDroppedP4[i] = candidate->SoftDroppedP4[i];
[de6d698]851 }
[9040259]852
[e9c0d73]853 //--- exclusive clustering variables ---
854 entry->ExclYmerge23 = candidate->ExclYmerge23;
855 entry->ExclYmerge34 = candidate->ExclYmerge34;
856 entry->ExclYmerge45 = candidate->ExclYmerge45;
[341014c]857 entry->ExclYmerge56 = candidate->ExclYmerge56;
[e9c0d73]858
[d7d2da3]859 FillParticles(candidate, &entry->Particles);
860 }
861}
862
863//------------------------------------------------------------------------------
864
865void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
866{
867 Candidate *candidate = 0;
868 MissingET *entry = 0;
869
870 // get the first entry
[341014c]871 if((candidate = static_cast<Candidate *>(array->At(0))))
[d7d2da3]872 {
873 const TLorentzVector &momentum = candidate->Momentum;
874
[341014c]875 entry = static_cast<MissingET *>(branch->NewEntry());
[d7d2da3]876
[3b465ca]877 entry->Eta = (-momentum).Eta();
[d7d2da3]878 entry->Phi = (-momentum).Phi();
879 entry->MET = momentum.Pt();
880 }
881}
882
883//------------------------------------------------------------------------------
884
885void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
886{
887 Candidate *candidate = 0;
888 ScalarHT *entry = 0;
889
890 // get the first entry
[341014c]891 if((candidate = static_cast<Candidate *>(array->At(0))))
[d7d2da3]892 {
893 const TLorentzVector &momentum = candidate->Momentum;
894
[341014c]895 entry = static_cast<ScalarHT *>(branch->NewEntry());
[d7d2da3]896
897 entry->HT = momentum.Pt();
898 }
899}
900
901//------------------------------------------------------------------------------
902
[71648c2]903void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
904{
[3b465ca]905 TIter iterator(array);
[71648c2]906 Candidate *candidate = 0;
907 Rho *entry = 0;
908
[3b465ca]909 // loop over all rho
910 iterator.Reset();
[341014c]911 while((candidate = static_cast<Candidate *>(iterator.Next())))
[71648c2]912 {
913 const TLorentzVector &momentum = candidate->Momentum;
914
[341014c]915 entry = static_cast<Rho *>(branch->NewEntry());
[71648c2]916
[8608ef8]917 entry->Rho = momentum.E();
[3b465ca]918 entry->Edges[0] = candidate->Edges[0];
919 entry->Edges[1] = candidate->Edges[1];
[71648c2]920 }
921}
922
923//------------------------------------------------------------------------------
924
[2e229c9]925void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
926{
927 Candidate *candidate = 0;
928 Weight *entry = 0;
929
930 // get the first entry
[341014c]931 if((candidate = static_cast<Candidate *>(array->At(0))))
[2e229c9]932 {
933 const TLorentzVector &momentum = candidate->Momentum;
934
[341014c]935 entry = static_cast<Weight *>(branch->NewEntry());
[2e229c9]936
937 entry->Weight = momentum.E();
938 }
939}
940
941//------------------------------------------------------------------------------
942
[8f7db23]943void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
944{
945 TIter iterator(array);
946 Candidate *candidate = 0;
947 HectorHit *entry = 0;
948
949 // loop over all roman pot hits
950 iterator.Reset();
[341014c]951 while((candidate = static_cast<Candidate *>(iterator.Next())))
[8f7db23]952 {
953 const TLorentzVector &position = candidate->Position;
954 const TLorentzVector &momentum = candidate->Momentum;
955
[341014c]956 entry = static_cast<HectorHit *>(branch->NewEntry());
[8f7db23]957
958 entry->E = momentum.E();
959
960 entry->Tx = momentum.Px();
961 entry->Ty = momentum.Py();
962
963 entry->T = position.T();
964
965 entry->X = position.X();
966 entry->Y = position.Y();
967 entry->S = position.Z();
[64a4950]968
969 entry->Particle = candidate->GetCandidates()->At(0);
[8f7db23]970 }
971}
972
973//------------------------------------------------------------------------------
974
[d7d2da3]975void TreeWriter::Process()
976{
977 TBranchMap::iterator itBranchMap;
978 ExRootTreeBranch *branch;
979 TProcessMethod method;
980 TObjArray *array;
981
982 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
983 {
984 branch = itBranchMap->first;
985 method = itBranchMap->second.first;
986 array = itBranchMap->second.second;
987
988 (this->*method)(branch, array);
989 }
990}
991
992//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.