Fork me on GitHub

source: git/modules/TreeWriter.cc@ 914d1be

Last change on this file since 914d1be was 781af69, checked in by michele <michele.selvaggi@…>, 3 years ago

added dNdx member

  • Property mode set to 100644
File size: 29.0 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;
[781af69]358 entry->dNdx = candidate->dNdx;
[2a3eb22]359
[341014c]360 entry->ErrorP = candidate->ErrorP;
361 entry->ErrorPT = candidate->ErrorPT;
[2671df6]362
363 // diagonal covariance matrix terms
364 entry->ErrorD0 = candidate->ErrorD0;
365 entry->ErrorC = candidate->ErrorC;
[341014c]366 entry->ErrorPhi = candidate->ErrorPhi;
[2671df6]367 entry->ErrorDZ = candidate->ErrorDZ;
368 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
369
370 // add some offdiagonal covariance matrix elements
371 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
372 entry->ErrorD0C = candidate->TrackCovariance(0,2);
373 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
374 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
375 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
376 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
377 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
378 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
379 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
380 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
[5496767]381
[e4c3fef]382 entry->Xd = candidate->Xd;
383 entry->Yd = candidate->Yd;
384 entry->Zd = candidate->Zd;
[9040259]385
[d7d2da3]386 const TLorentzVector &momentum = candidate->Momentum;
387
388 pt = momentum.Pt();
[2a3eb22]389 p = momentum.P();
390 phi = momentum.Phi();
[fd4b326]391 m = momentum.M();
[341014c]392 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
[2a3eb22]393
[d7d2da3]394 cosTheta = TMath::Abs(momentum.CosTheta());
395 signz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]396 eta = (cosTheta == 1.0 ? signz * 999.9 : momentum.Eta());
397 rapidity = (cosTheta == 1.0 ? signz * 999.9 : momentum.Rapidity());
[d7d2da3]398
[c1ea422]399 entry->P = p;
[341014c]400 entry->PT = pt;
[d7d2da3]401 entry->Eta = eta;
[2a3eb22]402 entry->Phi = phi;
403 entry->CtgTheta = ctgTheta;
[17cd992]404 entry->C = candidate->C;
[fd4b326]405 entry->Mass = m;
[5496767]406
[0060caa]407 //particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
408 //const TLorentzVector &initialPosition = particle->Position;
409 const TLorentzVector &initialPosition = candidate->InitialPosition;
[d7d2da3]410
411 entry->X = initialPosition.X();
412 entry->Y = initialPosition.Y();
413 entry->Z = initialPosition.Z();
[341014c]414 entry->T = initialPosition.T() * 1.0E-3 / c_light;
[d7d2da3]415
416 entry->Particle = particle;
[2600216]417
418 entry->VertexIndex = candidate->ClusterIndex;
[d7d2da3]419 }
420}
421
422//------------------------------------------------------------------------------
423
424void TreeWriter::ProcessTowers(ExRootTreeBranch *branch, TObjArray *array)
425{
426 TIter iterator(array);
427 Candidate *candidate = 0;
428 Tower *entry = 0;
429 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]430 const Double_t c_light = 2.99792458E8;
[8f7db23]431
[d07e957]432 // loop over all towers
[d7d2da3]433 iterator.Reset();
[341014c]434 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]435 {
436 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]437 const TLorentzVector &position = candidate->Position;
[8f7db23]438
[d7d2da3]439 pt = momentum.Pt();
440 cosTheta = TMath::Abs(momentum.CosTheta());
441 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]442 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
443 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]444
[341014c]445 entry = static_cast<Tower *>(branch->NewEntry());
[d7d2da3]446
447 entry->SetBit(kIsReferenced);
448 entry->SetUniqueID(candidate->GetUniqueID());
449
450 entry->Eta = eta;
451 entry->Phi = momentum.Phi();
452 entry->ET = pt;
453 entry->E = momentum.E();
454 entry->Eem = candidate->Eem;
455 entry->Ehad = candidate->Ehad;
456 entry->Edges[0] = candidate->Edges[0];
457 entry->Edges[1] = candidate->Edges[1];
458 entry->Edges[2] = candidate->Edges[2];
459 entry->Edges[3] = candidate->Edges[3];
[8f7db23]460
[341014c]461 entry->T = position.T() * 1.0E-3 / c_light;
[839deb7]462 entry->NTimeHits = candidate->NTimeHits;
[8f7db23]463
[d7d2da3]464 FillParticles(candidate, &entry->Particles);
465 }
466}
467
468//------------------------------------------------------------------------------
469
[4d7014e]470void TreeWriter::ProcessParticleFlowCandidates(ExRootTreeBranch *branch, TObjArray *array)
471{
472
473 TIter iterator(array);
474 Candidate *candidate = 0;
475 Candidate *particle = 0;
476 ParticleFlowCandidate *entry = 0;
[fd4b326]477 Double_t e, pt, signz, cosTheta, eta, rapidity, p, ctgTheta, phi, m;
[4d7014e]478 const Double_t c_light = 2.99792458E8;
479
480 // loop over all tracks
481 iterator.Reset();
482 while((candidate = static_cast<Candidate *>(iterator.Next())))
483 {
484 const TLorentzVector &position = candidate->Position;
485
486 cosTheta = TMath::Abs(position.CosTheta());
487 signz = (position.Pz() >= 0.0) ? 1.0 : -1.0;
488 eta = (cosTheta == 1.0 ? signz * 999.9 : position.Eta());
489 rapidity = (cosTheta == 1.0 ? signz * 999.9 : position.Rapidity());
490
491 entry = static_cast<ParticleFlowCandidate *>(branch->NewEntry());
492
493 entry->SetBit(kIsReferenced);
494 entry->SetUniqueID(candidate->GetUniqueID());
495
496 entry->PID = candidate->PID;
497
498 entry->Charge = candidate->Charge;
499
500 entry->EtaOuter = eta;
501 entry->PhiOuter = position.Phi();
502
503 entry->XOuter = position.X();
504 entry->YOuter = position.Y();
505 entry->ZOuter = position.Z();
506 entry->TOuter = position.T() * 1.0E-3 / c_light;
507
508 entry->L = candidate->L;
509
510 entry->D0 = candidate->D0;
511 entry->DZ = candidate->DZ;
[a95da74]512 entry->Nclusters = candidate->Nclusters;
[781af69]513 entry->dNdx = candidate->dNdx;
[4d7014e]514
515 entry->ErrorP = candidate->ErrorP;
516 entry->ErrorPT = candidate->ErrorPT;
517 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
[2671df6]518
519
520 // diagonal covariance matrix terms
521
522 entry->ErrorD0 = candidate->ErrorD0;
523 entry->ErrorC = candidate->ErrorC;
[4d7014e]524 entry->ErrorPhi = candidate->ErrorPhi;
[2671df6]525 entry->ErrorDZ = candidate->ErrorDZ;
526 entry->ErrorCtgTheta = candidate->ErrorCtgTheta;
527
528 // add some offdiagonal covariance matrix elements
529 entry->ErrorD0Phi = candidate->TrackCovariance(0,1);
530 entry->ErrorD0C = candidate->TrackCovariance(0,2);
531 entry->ErrorD0DZ = candidate->TrackCovariance(0,3);
532 entry->ErrorD0CtgTheta = candidate->TrackCovariance(0,4);
533 entry->ErrorPhiC = candidate->TrackCovariance(1,2);
534 entry->ErrorPhiDZ = candidate->TrackCovariance(1,3);
535 entry->ErrorPhiCtgTheta = candidate->TrackCovariance(1,4);
536 entry->ErrorCDZ = candidate->TrackCovariance(2,3);
537 entry->ErrorCCtgTheta = candidate->TrackCovariance(2,4);
538 entry->ErrorDZCtgTheta = candidate->TrackCovariance(3,4);
539
[4d7014e]540 entry->Xd = candidate->Xd;
541 entry->Yd = candidate->Yd;
542 entry->Zd = candidate->Zd;
543
544 const TLorentzVector &momentum = candidate->Momentum;
545
546 e = momentum.E();
547 pt = momentum.Pt();
548 p = momentum.P();
549 phi = momentum.Phi();
[fd4b326]550 m = momentum.M();
[4d7014e]551 ctgTheta = (TMath::Tan(momentum.Theta()) != 0) ? 1 / TMath::Tan(momentum.Theta()) : 1e10;
552
553 entry->E = e;
554 entry->P = p;
555 entry->PT = pt;
556 entry->Eta = eta;
557 entry->Phi = phi;
558 entry->CtgTheta = ctgTheta;
[17cd992]559 entry->C = candidate->C;
[fd4b326]560 entry->Mass = m;
[dc883b4]561
[0060caa]562 //particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
563 //const TLorentzVector &initialPosition = particle->Position;
564 const TLorentzVector &initialPosition = candidate->InitialPosition;
[4d7014e]565
566 entry->X = initialPosition.X();
567 entry->Y = initialPosition.Y();
568 entry->Z = initialPosition.Z();
569 entry->T = initialPosition.T() * 1.0E-3 / c_light;
570
571 entry->VertexIndex = candidate->ClusterIndex;
572
573 entry->Eem = candidate->Eem;
574 entry->Ehad = candidate->Ehad;
575 entry->Edges[0] = candidate->Edges[0];
576 entry->Edges[1] = candidate->Edges[1];
577 entry->Edges[2] = candidate->Edges[2];
578 entry->Edges[3] = candidate->Edges[3];
579
[0060caa]580 //entry->T = position.T() * 1.0E-3 / c_light;
[4d7014e]581 entry->NTimeHits = candidate->NTimeHits;
582
583 FillParticles(candidate, &entry->Particles);
584
585 }
586}
587
588//------------------------------------------------------------------------------
589
[d7d2da3]590void TreeWriter::ProcessPhotons(ExRootTreeBranch *branch, TObjArray *array)
591{
592 TIter iterator(array);
593 Candidate *candidate = 0;
594 Photon *entry = 0;
595 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]596 const Double_t c_light = 2.99792458E8;
[8f7db23]597
[d7d2da3]598 array->Sort();
599
600 // loop over all photons
601 iterator.Reset();
[341014c]602 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]603 {
604 TIter it1(candidate->GetCandidates());
605 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]606 const TLorentzVector &position = candidate->Position;
[8f7db23]607
[d7d2da3]608 pt = momentum.Pt();
609 cosTheta = TMath::Abs(momentum.CosTheta());
610 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]611 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
612 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]613
[341014c]614 entry = static_cast<Photon *>(branch->NewEntry());
[d7d2da3]615
616 entry->Eta = eta;
617 entry->Phi = momentum.Phi();
618 entry->PT = pt;
619 entry->E = momentum.E();
[341014c]620 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]621
[d5cae2b]622 // Isolation variables
[9040259]623
[d5cae2b]624 entry->IsolationVar = candidate->IsolationVar;
[341014c]625 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
626 entry->SumPtCharged = candidate->SumPtCharged;
627 entry->SumPtNeutral = candidate->SumPtNeutral;
628 entry->SumPtChargedPU = candidate->SumPtChargedPU;
629 entry->SumPt = candidate->SumPt;
[d5cae2b]630
[341014c]631 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad / candidate->Eem : 999.9;
[d7d2da3]632
[0e0f211]633 // 1: prompt -- 2: non prompt -- 3: fake
634 entry->Status = candidate->Status;
635
[d7d2da3]636 FillParticles(candidate, &entry->Particles);
637 }
638}
639
640//------------------------------------------------------------------------------
641
642void TreeWriter::ProcessElectrons(ExRootTreeBranch *branch, TObjArray *array)
643{
644 TIter iterator(array);
645 Candidate *candidate = 0;
646 Electron *entry = 0;
647 Double_t pt, signPz, cosTheta, eta, rapidity;
[22dc7fd]648 const Double_t c_light = 2.99792458E8;
[8f7db23]649
[d7d2da3]650 array->Sort();
651
652 // loop over all electrons
653 iterator.Reset();
[341014c]654 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]655 {
656 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]657 const TLorentzVector &position = candidate->Position;
[8f7db23]658
[d7d2da3]659 pt = momentum.Pt();
660 cosTheta = TMath::Abs(momentum.CosTheta());
661 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]662 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
663 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]664
[341014c]665 entry = static_cast<Electron *>(branch->NewEntry());
[d7d2da3]666
667 entry->Eta = eta;
668 entry->Phi = momentum.Phi();
669 entry->PT = pt;
[8f7db23]670
[341014c]671 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]672
[0518688]673 // displacement
[341014c]674 entry->D0 = candidate->D0;
675 entry->ErrorD0 = candidate->ErrorD0;
676 entry->DZ = candidate->DZ;
677 entry->ErrorDZ = candidate->ErrorDZ;
[9040259]678
[0518688]679 // Isolation variables
[d5cae2b]680 entry->IsolationVar = candidate->IsolationVar;
[341014c]681 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
682 entry->SumPtCharged = candidate->SumPtCharged;
683 entry->SumPtNeutral = candidate->SumPtNeutral;
684 entry->SumPtChargedPU = candidate->SumPtChargedPU;
685 entry->SumPt = candidate->SumPt;
[d5cae2b]686
[d7d2da3]687 entry->Charge = candidate->Charge;
688
689 entry->EhadOverEem = 0.0;
690
691 entry->Particle = candidate->GetCandidates()->At(0);
692 }
693}
694
695//------------------------------------------------------------------------------
696
697void TreeWriter::ProcessMuons(ExRootTreeBranch *branch, TObjArray *array)
698{
699 TIter iterator(array);
700 Candidate *candidate = 0;
701 Muon *entry = 0;
702 Double_t pt, signPz, cosTheta, eta, rapidity;
[8f7db23]703
[22dc7fd]704 const Double_t c_light = 2.99792458E8;
[8f7db23]705
[d7d2da3]706 array->Sort();
707
708 // loop over all muons
709 iterator.Reset();
[341014c]710 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]711 {
712 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]713 const TLorentzVector &position = candidate->Position;
[8f7db23]714
[d7d2da3]715 pt = momentum.Pt();
716 cosTheta = TMath::Abs(momentum.CosTheta());
717 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]718 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
719 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]720
[341014c]721 entry = static_cast<Muon *>(branch->NewEntry());
[d7d2da3]722
723 entry->SetBit(kIsReferenced);
724 entry->SetUniqueID(candidate->GetUniqueID());
725
726 entry->Eta = eta;
727 entry->Phi = momentum.Phi();
728 entry->PT = pt;
729
[341014c]730 entry->T = position.T() * 1.0E-3 / c_light;
[0518688]731
732 // displacement
[341014c]733 entry->D0 = candidate->D0;
734 entry->ErrorD0 = candidate->ErrorD0;
735 entry->DZ = candidate->DZ;
736 entry->ErrorDZ = candidate->ErrorDZ;
[0518688]737
[d5cae2b]738 // Isolation variables
[9040259]739
[d5cae2b]740 entry->IsolationVar = candidate->IsolationVar;
[341014c]741 entry->IsolationVarRhoCorr = candidate->IsolationVarRhoCorr;
742 entry->SumPtCharged = candidate->SumPtCharged;
743 entry->SumPtNeutral = candidate->SumPtNeutral;
744 entry->SumPtChargedPU = candidate->SumPtChargedPU;
745 entry->SumPt = candidate->SumPt;
[8f7db23]746
[d7d2da3]747 entry->Charge = candidate->Charge;
748
749 entry->Particle = candidate->GetCandidates()->At(0);
750 }
751}
752
753//------------------------------------------------------------------------------
754
755void TreeWriter::ProcessJets(ExRootTreeBranch *branch, TObjArray *array)
756{
757 TIter iterator(array);
758 Candidate *candidate = 0, *constituent = 0;
759 Jet *entry = 0;
760 Double_t pt, signPz, cosTheta, eta, rapidity;
761 Double_t ecalEnergy, hcalEnergy;
[22dc7fd]762 const Double_t c_light = 2.99792458E8;
[de6d698]763 Int_t i;
[8f7db23]764
[d7d2da3]765 array->Sort();
766
767 // loop over all jets
768 iterator.Reset();
[341014c]769 while((candidate = static_cast<Candidate *>(iterator.Next())))
[d7d2da3]770 {
771 TIter itConstituents(candidate->GetCandidates());
[8f7db23]772
[d7d2da3]773 const TLorentzVector &momentum = candidate->Momentum;
[22dc7fd]774 const TLorentzVector &position = candidate->Position;
[8f7db23]775
[d7d2da3]776 pt = momentum.Pt();
777 cosTheta = TMath::Abs(momentum.CosTheta());
778 signPz = (momentum.Pz() >= 0.0) ? 1.0 : -1.0;
[341014c]779 eta = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Eta());
780 rapidity = (cosTheta == 1.0 ? signPz * 999.9 : momentum.Rapidity());
[d7d2da3]781
[341014c]782 entry = static_cast<Jet *>(branch->NewEntry());
[d7d2da3]783
784 entry->Eta = eta;
785 entry->Phi = momentum.Phi();
786 entry->PT = pt;
787
[341014c]788 entry->T = position.T() * 1.0E-3 / c_light;
[8f7db23]789
[d7d2da3]790 entry->Mass = momentum.M();
791
[ba1f1ee]792 entry->Area = candidate->Area;
793
[d7d2da3]794 entry->DeltaEta = candidate->DeltaEta;
795 entry->DeltaPhi = candidate->DeltaPhi;
796
[fe0273c]797 entry->Flavor = candidate->Flavor;
798 entry->FlavorAlgo = candidate->FlavorAlgo;
799 entry->FlavorPhys = candidate->FlavorPhys;
800
[d7d2da3]801 entry->BTag = candidate->BTag;
[9040259]802
803 entry->BTagAlgo = candidate->BTagAlgo;
[fe0273c]804 entry->BTagPhys = candidate->BTagPhys;
[9040259]805
[d7d2da3]806 entry->TauTag = candidate->TauTag;
[7429c6a]807 entry->TauWeight = candidate->TauWeight;
[d7d2da3]808
809 entry->Charge = candidate->Charge;
810
811 itConstituents.Reset();
812 entry->Constituents.Clear();
813 ecalEnergy = 0.0;
814 hcalEnergy = 0.0;
[341014c]815 while((constituent = static_cast<Candidate *>(itConstituents.Next())))
[d7d2da3]816 {
817 entry->Constituents.Add(constituent);
818 ecalEnergy += constituent->Eem;
819 hcalEnergy += constituent->Ehad;
820 }
821
[341014c]822 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy / ecalEnergy : 999.9;
[8f7db23]823
[24d005f]824 //--- Pile-Up Jet ID variables ----
825
826 entry->NCharged = candidate->NCharged;
827 entry->NNeutrals = candidate->NNeutrals;
[c614dd7]828
829 entry->NeutralEnergyFraction = candidate->NeutralEnergyFraction;
830 entry->ChargedEnergyFraction = candidate->ChargedEnergyFraction;
[24d005f]831 entry->Beta = candidate->Beta;
832 entry->BetaStar = candidate->BetaStar;
833 entry->MeanSqDeltaR = candidate->MeanSqDeltaR;
834 entry->PTD = candidate->PTD;
[9040259]835
[de6d698]836 //--- Sub-structure variables ----
[9040259]837
838 entry->NSubJetsTrimmed = candidate->NSubJetsTrimmed;
839 entry->NSubJetsPruned = candidate->NSubJetsPruned;
[de6d698]840 entry->NSubJetsSoftDropped = candidate->NSubJetsSoftDropped;
[9040259]841
[341014c]842 entry->SoftDroppedJet = candidate->SoftDroppedJet;
843 entry->SoftDroppedSubJet1 = candidate->SoftDroppedSubJet1;
[ba75867]844 entry->SoftDroppedSubJet2 = candidate->SoftDroppedSubJet2;
845
[de6d698]846 for(i = 0; i < 5; i++)
847 {
[341014c]848 entry->FracPt[i] = candidate->FracPt[i];
849 entry->Tau[i] = candidate->Tau[i];
850 entry->TrimmedP4[i] = candidate->TrimmedP4[i];
851 entry->PrunedP4[i] = candidate->PrunedP4[i];
852 entry->SoftDroppedP4[i] = candidate->SoftDroppedP4[i];
[de6d698]853 }
[9040259]854
[e9c0d73]855 //--- exclusive clustering variables ---
856 entry->ExclYmerge23 = candidate->ExclYmerge23;
857 entry->ExclYmerge34 = candidate->ExclYmerge34;
858 entry->ExclYmerge45 = candidate->ExclYmerge45;
[341014c]859 entry->ExclYmerge56 = candidate->ExclYmerge56;
[e9c0d73]860
[d7d2da3]861 FillParticles(candidate, &entry->Particles);
862 }
863}
864
865//------------------------------------------------------------------------------
866
867void TreeWriter::ProcessMissingET(ExRootTreeBranch *branch, TObjArray *array)
868{
869 Candidate *candidate = 0;
870 MissingET *entry = 0;
871
872 // get the first entry
[341014c]873 if((candidate = static_cast<Candidate *>(array->At(0))))
[d7d2da3]874 {
875 const TLorentzVector &momentum = candidate->Momentum;
876
[341014c]877 entry = static_cast<MissingET *>(branch->NewEntry());
[d7d2da3]878
[3b465ca]879 entry->Eta = (-momentum).Eta();
[d7d2da3]880 entry->Phi = (-momentum).Phi();
881 entry->MET = momentum.Pt();
882 }
883}
884
885//------------------------------------------------------------------------------
886
887void TreeWriter::ProcessScalarHT(ExRootTreeBranch *branch, TObjArray *array)
888{
889 Candidate *candidate = 0;
890 ScalarHT *entry = 0;
891
892 // get the first entry
[341014c]893 if((candidate = static_cast<Candidate *>(array->At(0))))
[d7d2da3]894 {
895 const TLorentzVector &momentum = candidate->Momentum;
896
[341014c]897 entry = static_cast<ScalarHT *>(branch->NewEntry());
[d7d2da3]898
899 entry->HT = momentum.Pt();
900 }
901}
902
903//------------------------------------------------------------------------------
904
[71648c2]905void TreeWriter::ProcessRho(ExRootTreeBranch *branch, TObjArray *array)
906{
[3b465ca]907 TIter iterator(array);
[71648c2]908 Candidate *candidate = 0;
909 Rho *entry = 0;
910
[3b465ca]911 // loop over all rho
912 iterator.Reset();
[341014c]913 while((candidate = static_cast<Candidate *>(iterator.Next())))
[71648c2]914 {
915 const TLorentzVector &momentum = candidate->Momentum;
916
[341014c]917 entry = static_cast<Rho *>(branch->NewEntry());
[71648c2]918
[8608ef8]919 entry->Rho = momentum.E();
[3b465ca]920 entry->Edges[0] = candidate->Edges[0];
921 entry->Edges[1] = candidate->Edges[1];
[71648c2]922 }
923}
924
925//------------------------------------------------------------------------------
926
[2e229c9]927void TreeWriter::ProcessWeight(ExRootTreeBranch *branch, TObjArray *array)
928{
929 Candidate *candidate = 0;
930 Weight *entry = 0;
931
932 // get the first entry
[341014c]933 if((candidate = static_cast<Candidate *>(array->At(0))))
[2e229c9]934 {
935 const TLorentzVector &momentum = candidate->Momentum;
936
[341014c]937 entry = static_cast<Weight *>(branch->NewEntry());
[2e229c9]938
939 entry->Weight = momentum.E();
940 }
941}
942
943//------------------------------------------------------------------------------
944
[8f7db23]945void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array)
946{
947 TIter iterator(array);
948 Candidate *candidate = 0;
949 HectorHit *entry = 0;
950
951 // loop over all roman pot hits
952 iterator.Reset();
[341014c]953 while((candidate = static_cast<Candidate *>(iterator.Next())))
[8f7db23]954 {
955 const TLorentzVector &position = candidate->Position;
956 const TLorentzVector &momentum = candidate->Momentum;
957
[341014c]958 entry = static_cast<HectorHit *>(branch->NewEntry());
[8f7db23]959
960 entry->E = momentum.E();
961
962 entry->Tx = momentum.Px();
963 entry->Ty = momentum.Py();
964
965 entry->T = position.T();
966
967 entry->X = position.X();
968 entry->Y = position.Y();
969 entry->S = position.Z();
[64a4950]970
971 entry->Particle = candidate->GetCandidates()->At(0);
[8f7db23]972 }
973}
974
975//------------------------------------------------------------------------------
976
[d7d2da3]977void TreeWriter::Process()
978{
979 TBranchMap::iterator itBranchMap;
980 ExRootTreeBranch *branch;
981 TProcessMethod method;
982 TObjArray *array;
983
984 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
985 {
986 branch = itBranchMap->first;
987 method = itBranchMap->second.first;
988 array = itBranchMap->second.second;
989
990 (this->*method)(branch, array);
991 }
992}
993
994//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.