ATLAS Offline Software
Loading...
Searching...
No Matches
PhysicsAnalysis/TauID/TauAnalysisTools/Root/HelperFunctions.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <fstream>
6
7// local include(s)
11#include "TF1.h"
12
13#ifdef ASGTOOL_ATHENA
14#include "CLHEP/Units/SystemOfUnits.h"
15using CLHEP::GeV;
16#else
17#define GeV 1000
18#endif
19
20using namespace TauAnalysisTools;
21
22//______________________________________________________________________________
23void TauAnalysisTools::split(const std::string& sInput, const char cDelim, std::vector<std::string>& vOut)
24{
25 std::stringstream sSS(sInput);
26 std::string sItem;
27 while (std::getline(sSS, sItem, cDelim))
28 vOut.push_back(sItem);
29}
30
31//______________________________________________________________________________
32void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<std::string>& vOut)
33{
34 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
35 std::string sItem;
36 while (std::getline(sSS, sItem, cDelim))
37 vOut.push_back(sItem);
38}
39
40//______________________________________________________________________________
41void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<size_t>& vOut)
42{
43 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
44 std::string sItem;
45 while (std::getline(sSS, sItem, cDelim))
46 vOut.push_back(size_t(stoi(sItem)));
47}
48
49//______________________________________________________________________________
50void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<int>& vOut)
51{
52 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
53 std::string sItem;
54 while (std::getline(sSS, sItem, cDelim))
55 vOut.push_back(stoi(sItem));
56}
57
58//______________________________________________________________________________
59void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<unsigned>& vOut)
60{
61 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
62 std::string sItem;
63 while (std::getline(sSS, sItem, cDelim))
64 vOut.push_back(stoi(sItem));
65}
66
67//______________________________________________________________________________
68void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<float>& vOut)
69{
70 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
71 std::string sItem;
72 while (std::getline(sSS, sItem, cDelim))
73 vOut.push_back(stof(sItem));
74}
75
76//______________________________________________________________________________
77void TauAnalysisTools::split(TEnv& rEnv, const std::string& sKey, const char cDelim, std::vector<double>& vOut)
78{
79 std::stringstream sSS(rEnv.GetValue(sKey.c_str(),""));
80 std::string sItem;
81 while (std::getline(sSS, sItem, cDelim))
82 vOut.push_back(stod(sItem));
83}
84//______________________________________________________________________________
86{
87 // return tau pt in GeV
88 return xTau.pt()/1000.;
89}
90
91//______________________________________________________________________________
93{
94 // return tau P in GeV
95 return xTau.p4().P()/1000.;
96}
97
98//______________________________________________________________________________
100{
101 // return tau eta
102 return xTau.eta();
103}
104
105//______________________________________________________________________________
107{
108 // return absolute tau eta
109 return std::abs(xTau.eta());
110}
111
112//______________________________________________________________________________
114{
115 // return MVA based tau pt in GeV
116 return xTau.ptFinalCalib()/GeV;
117}
118
119//______________________________________________________________________________
121{
122 // return MVA based tau eta
123 return xTau.etaFinalCalib();
124}
125
126//______________________________________________________________________________
128{
129 // return MVA based absolute tau eta
130 return std::abs(xTau.etaFinalCalib());
131}
132
133//______________________________________________________________________________
135{
136 // return tau P in GeV
138}
139
140//______________________________________________________________________________
142{
143 // return leading charge tau track eta
144 double dTrackEta = 0.;
145 double dTrackMaxPt = 0.;
146 for( unsigned int iNumTrack = 0; iNumTrack < xTau.nTracks(); iNumTrack++)
147 {
148 if (xTau.track(iNumTrack)->pt() > dTrackMaxPt)
149 {
150 dTrackMaxPt = xTau.track(iNumTrack)->pt();
151 dTrackEta = xTau.track(iNumTrack)->eta();
152 }
153 }
154 return dTrackEta;
155}
156
157//______________________________________________________________________________
159{
160 // return truth tau Pt in GeV
161 const xAOD::TruthParticle* xTruthTau = getTruth(xTau);
162
163 // if there is a truth tau return pT, otherwise return 0 (getTruth will print an error)
164 static const SG::ConstAccessor<char> accIsHadronicTau ("IsHadronicTau");
165 if (xTruthTau!=nullptr && accIsHadronicTau (*xTruthTau))
166 return xTruthTau->pt()/GeV;
167 else
168 return 0.;
169}
170//______________________________________________________________________________
172{
173 // return truth visible tau Pt in GeV
174 const xAOD::TruthParticle* xTruthTau = getTruth(xTau);
175
176 // if there is a truth tau return visible pT, otherwise return 0 (getTruth will print an error)
177 static const SG::ConstAccessor<char> accIsHadronicTau ("IsHadronicTau");
178 static const SG::ConstAccessor<double> accPtVis("pt_vis");
179 if (xTruthTau!=nullptr && accIsHadronicTau (*xTruthTau))
180 return (accPtVis(*xTruthTau)/GeV);
181 else
182 return 0.;
183}
184//______________________________________________________________________________
186{
187 // return truth tau absolute eta
188 const xAOD::TruthParticle* xTruthTau = getTruth(xTau);
189
190 // if there is a truth tau return absolute eta, otherwise return -5 (getTruth will print an error)
191 static const SG::ConstAccessor<char> accIsHadronicTau ("IsHadronicTau");
192 if (xTruthTau!=nullptr && accIsHadronicTau (*xTruthTau))
193 return std::abs(xTruthTau->eta());
194 else
195 return -5.;
196}
197
198//______________________________________________________________________________
200{
201 // return truth tau decay mode.
202 int iDecayMode = getTruthDecayMode(xTau);
203 return static_cast<double>(iDecayMode);
204}
205
206//______________________________________________________________________________
208{
210 static const SG::ConstAccessor<Link_t> accTruthParticleLink("truthParticleLink");
211 if (!accTruthParticleLink(xTau))
212 {
213 Error("TauAnalysisTools::getTruth", "No truth match information available. Please run TauTruthMatchingTool first");
214 }
215
216 const Link_t xTruthTauLink = accTruthParticleLink(xTau);
217 const xAOD::TruthParticle* xTruthTau = xTruthTauLink.cachedElement();
218
219 return xTruthTau;
220}
221
222
223//______________________________________________________________________________
225{
226 const xAOD::TruthParticle* xTruthTau = getTruth(xTau);
227
228 static const SG::ConstAccessor<char> accIsHadronicTau ("IsHadronicTau");
229 if (xTruthTau!=nullptr && accIsHadronicTau(*xTruthTau))
230 return getTruthDecayMode(*xTruthTau);
231 else
233}
234
235//______________________________________________________________________________
237{
238 static const SG::ConstAccessor<size_t> accNumCharged ("numCharged");
239 if (!(accNumCharged.isAvailable(xTruthTau)))
240 {
241 // passed truth particle is not a truth tau
243 }
244
245 int iCharged = getNTauDecayParticles(xTruthTau,MC::PIPLUS, true) + getNTauDecayParticles(xTruthTau,MC::KPLUS, true);
246 int iNeutral = getNTauDecayParticles(xTruthTau,MC::PI0, true);
247 if (iCharged == 1)
248 {
249 if (iNeutral == 0) return xAOD::TauJetParameters::DecayMode::Mode_1p0n;
250 if (iNeutral == 1) return xAOD::TauJetParameters::DecayMode::Mode_1p1n;
251 if (iNeutral >= 2) return xAOD::TauJetParameters::DecayMode::Mode_1pXn;
252 }
253 else if (iCharged == 3)
254 {
255 if (iNeutral == 0) return xAOD::TauJetParameters::DecayMode::Mode_3p0n;
256 if (iNeutral >= 1) return xAOD::TauJetParameters::DecayMode::Mode_3pXn;
257 }
258
259 if (iCharged == 2 or iCharged == 4 or iCharged == 5)
261 if (iCharged == 0 or iCharged >=6)
263
264 // if you got here, something should have gone wrong
266}
267
268//______________________________________________________________________________
269int TauAnalysisTools::getNTauDecayParticles(const xAOD::TruthParticle& xTruthTau, int iPdgId, bool bCompareAbsoluteValues)
270{
271 int iNum = 0;
272 static const SG::ConstAccessor<std::vector<int> > accDecayModeVector("DecayModeVector");
273 if (!accDecayModeVector.isAvailable(xTruthTau))
274 {
275 Warning("TauAnalysisTools::getNTauDecayParticles", "passed truth particle is not a truth tau, return 0");
276 return 0;
277 }
278
279 for(auto iPdgId2 : accDecayModeVector(xTruthTau))
280 if (!bCompareAbsoluteValues)
281 {
282 if (iPdgId2 == iPdgId) iNum++;
283 }
284 else
285 {
286 if (std::abs(iPdgId2) == std::abs(iPdgId)) iNum++;
287 }
288 return iNum;
289}
290
291//______________________________________________________________________________
293{
294 // returns true if last line in file is empty or the line starts with the
295 // number sign #
296
297 std::ifstream fInputFile;
298 fInputFile.open(sFileName);
299 if(!fInputFile.is_open())
300 return true;
301
302 fInputFile.seekg(-1,fInputFile.end);
303
304 bool bKeepLooping = true;
305 while(bKeepLooping)
306 {
307 char ch;
308 fInputFile.get(ch);
309
310 if(static_cast<int>(fInputFile.tellg()) <= 1)
311 {
312 fInputFile.seekg(0);
313 bKeepLooping = false;
314 }
315 else if(ch == '\n')
316 bKeepLooping = false;
317 else
318 fInputFile.seekg(-2,fInputFile.cur);
319 }
320
321 std::string sLastLine;
322 getline(fInputFile,sLastLine);
323 fInputFile.close();
324
325 return (sLastLine.size() == 0 or sLastLine[0] == '#');
326}
327
328//______________________________________________________________________________
329void TauAnalysisTools::createPi0Vectors(const xAOD::TauJet* xTau, std::vector<TLorentzVector>& vPi0s)
330{
331 // reset the pi0s
332 vPi0s.clear();
333
334 // Since the PFO links as they come out of reconstruction, only correspond to
335 // calorimeter clusters, whereas we want the consts_pi0 vectors to correspond
336 // to real pi0s, we need to be careful to collect the PFOs correctly to pi0s
337 // for the cases where number of pi0s does not match to the decay mode:
338 size_t iNumPi0PFO = xTau->nPi0PFOs();
339
340 int iDecayMode = -1;
341
343 {
344 Error("TauAnalysisTools::createPi0Vectors", "Failed to retrieve panTauDetail decay mode.");
345 return;
346 }
347
348 if (iDecayMode == xAOD::TauJetParameters::DecayMode::Mode_1p1n && iNumPi0PFO > 1)
349 {
350 // float fMassPi0 = ParticleConstants::piZeroMassInMeV;
352
353 // combine both photons (with 0 mass from Pantau) to one pi0 vector:
354 const xAOD::PFO* xPfo1 = xTau->pi0PFO(0);
355 const xAOD::PFO* xPfo2 = xTau->pi0PFO(1);
356 vPi0s.push_back(xPfo1->p4() + xPfo2->p4());
357
358 // re-set the mass to one pi0:
359 double dNewMomentum = std::sqrt(vPi0s[0].E() * vPi0s[0].E() - fMassPi0Squared);
360 vPi0s[0].SetPxPyPzE(vPi0s[0].Vect().Unit().Px() * dNewMomentum,
361 vPi0s[0].Vect().Unit().Py() * dNewMomentum,
362 vPi0s[0].Vect().Unit().Pz() * dNewMomentum,
363 vPi0s[0].E());
364 }
365 else if (iDecayMode == xAOD::TauJetParameters::DecayMode::Mode_1pXn && iNumPi0PFO == 1)
366 {
367 // make a single pi0 from a PFO that contains two pi0s:
368 const xAOD::PFO* xPfo = xTau->pi0PFO(0);
369 // add the 2-pi0 vector preliminarily to the pi0vector:
370 vPi0s.push_back(xPfo->p4());
371
372 // re-set the mass back to one pi0:
373 double dNewMomentum = std::sqrt(vPi0s[0].E() / 2 * vPi0s[0].E() / 2 - vPi0s[0].M() / 2. * vPi0s[0].M() / 2.);
374 vPi0s[0].SetVectM(vPi0s[0].Vect() * (dNewMomentum / vPi0s[0].P()), vPi0s[0].M() / 2.);
375
376 // create another pi0 from the same vector:
377 vPi0s.push_back(vPi0s[0]);
378 }
379 else
380 {
381 // if it's not any of the special cases above then just collect the PFOs:
382 for (size_t iPFO = 0; iPFO < iNumPi0PFO; iPFO++)
383 {
384 vPi0s.push_back(xTau->pi0PFO(iPFO)->p4());
385 }
386 }
387}
388
389
390//______________________________________________________________________________
391void TauAnalysisTools::correctedPi0Vectors(const xAOD::TauJet* xTau, std::vector<TLorentzVector>& correctedPi0s, TLorentzVector& TauP4){
392 //reset the pi0s
393 correctedPi0s.clear();
394
395 int iDecayMode = -1;
396
398 {
399 Error("TauAnalysisTools::correctedPi0Vectors", "Failed to retrieve panTauDetail decay mode.");
400 return;
401 }
402
403 //Reading in the pi0 vector from createPi0Vectors
404 std::vector<TLorentzVector> vPi0s;
405 createPi0Vectors(xTau,vPi0s);
406
408 //Adding up Pi0 P4s from createPi0Vectors
409 TLorentzVector Sum_vPi0s;
410 for(unsigned int i = 0; i < vPi0s.size() ; i++){
411 Sum_vPi0s += vPi0s[i];
412 }
413
414 //Get sum of the chargedPFO (i.e. tau track) p4
415 TLorentzVector Sum_ChrgPFOP4;
416 for(const xAOD::TauTrack* track : xTau->tracks()) {
417 Sum_ChrgPFOP4 += track->p4();
418 }
419
420 //Get tau FinalCalib P4 (explicitly requiring p4(xAOD::TauJetParameters::TauCalibType::FinalCalib) should be superfluous, as FinalCalib is the default p4)
421 TLorentzVector FinalCalibP4 = xTau->p4();
422
423 //Calculate the difference 3-vector between FinalCalib and Sum of chargedPFOP4
424 double px = FinalCalibP4.Px() - Sum_ChrgPFOP4.Px();
425 double py = FinalCalibP4.Py() - Sum_ChrgPFOP4.Py();
426 double pz = FinalCalibP4.Pz() - Sum_ChrgPFOP4.Pz();
427
428 double p_correctedPi0s = std::sqrt( std::pow(px,2.) + std::pow(py,2.) + std::pow(pz,2.) );
429 double p_vPi0s = Sum_vPi0s.P();
430
431 //Calucate scale factor for the pi0 3-vector momentum
432 double X = p_correctedPi0s/p_vPi0s;
433
434 //Scale the pi0s with X and recalculate the new pi0 energy
435 double px_scaled, py_scaled, pz_scaled, e;
437 for(unsigned int i = 0; i < vPi0s.size() ; i++){
438 px_scaled = vPi0s[i].Px() * X;
439 py_scaled = vPi0s[i].Py() * X;
440 pz_scaled = vPi0s[i].Pz() * X;
441 e = std::sqrt( std::pow(px_scaled,2.) + std::pow(py_scaled,2.) + std::pow(pz_scaled,2.) + std::pow(mPi0,2.) );
442
443 //Append the corrected pi0P4 to correctedPi0s
444 TLorentzVector P4_correctedPi0s;
445 P4_correctedPi0s.SetPxPyPzE(px_scaled,py_scaled,pz_scaled,e);
446 correctedPi0s.push_back(P4_correctedPi0s);
447 }
448 }else{
449 correctedPi0s = vPi0s;
450 }
451
452 //Correct angles between pi0s for 1pXn decays with 1 cluster
453 if(iDecayMode == xAOD::TauJetParameters::DecayMode::Mode_1pXn && xTau->nPi0PFOs() == 1){
454
455 //Get Function of Delta R between the two Pi0s
456 TF1 DeltaRdist("DeltaRdist", "pol3", 0, 67500);
457 DeltaRdist.SetParameter(0, 0.07924);
458 DeltaRdist.SetParameter(1, -2.078/1000000.);
459 DeltaRdist.SetParameter(2, 2.619/100000000000.);
460 DeltaRdist.SetParameter(3, -1.238/10000000000000000.);
461
462 //Get Sum of pi0 P4.Pt()
463 TLorentzVector SumPi0_P4;
464 for( unsigned int i = 0 ; i < correctedPi0s.size() ; i++){
465 SumPi0_P4 += correctedPi0s[i];
466 }
467
468 float SumPi0_pt = SumPi0_P4.Pt();
469
470 //Get delta R value (mean of true DeltaR distribution)
471 float deltaR;
472 if(SumPi0_pt >= 67500){
473 deltaR = 0.020; // = DeltaRdist.Eval(67500);
474 } else{
475 deltaR = DeltaRdist.Eval(SumPi0_pt);
476 }
477
478 TLorentzVector correctedPi0_0, correctedPi0_1;
479 correctedPi0_0.SetPtEtaPhiM( correctedPi0s[0].Pt()/cos(0.5*deltaR/std::sqrt(2.0)), correctedPi0s[0].Eta()+0.5*deltaR/std::sqrt(2.0), correctedPi0s[0].Phi()+0.5*deltaR/std::sqrt(2.0), correctedPi0s[0].M() );
480 correctedPi0_1.SetPtEtaPhiM( correctedPi0s[1].Pt()/cos(0.5*deltaR/std::sqrt(2.0)), correctedPi0s[1].Eta()-0.5*deltaR/std::sqrt(2.0), correctedPi0s[1].Phi()-0.5*deltaR/std::sqrt(2.0), correctedPi0s[1].M() );
481
482 std::vector<TLorentzVector> AngleCorrectedPi0s;
483 AngleCorrectedPi0s.push_back(correctedPi0_0);
484 AngleCorrectedPi0s.push_back(correctedPi0_1);
485
486 //Reparametrise: Delta R -> mass of pi0 Cluster
487 TLorentzVector PionCluster_angleCorrected = AngleCorrectedPi0s[0]+AngleCorrectedPi0s[1];
488
489 double dNewMomentum = std::sqrt(PionCluster_angleCorrected.E()/2. * PionCluster_angleCorrected.E()/2. - PionCluster_angleCorrected.M() / 2. * PionCluster_angleCorrected.M() / 2.);
490 correctedPi0s[0].SetVectM(PionCluster_angleCorrected.Vect() * (dNewMomentum / PionCluster_angleCorrected.P()), PionCluster_angleCorrected.M() / 2.);
491 correctedPi0s[1] = correctedPi0s[0];
492 }
493
494 //Calculate the new tau P4
495 for(const xAOD::TauTrack* track : xTau->tracks()) {
496 TauP4 += track->p4();
497 }
498
499 for(unsigned int iPi0=0; iPi0 < correctedPi0s.size(); iPi0++) {
500 TauP4 += correctedPi0s[iPi0];
501 }
502
503}
504
505//______________________________________________________________________________
507{
509 static const SG::ConstAccessor< Link_t > accTruthParticleLink("truthParticleLink");
510 if (!accTruthParticleLink.isAvailable(xTau))
511 Error("TauAnalysisTools::getTruthParticleType", "No truth match information available. Please run TauTruthMatchingTool first.");
512
513 const xAOD::TruthParticle* xTruthParticle = xAOD::TauHelpers::getTruthParticle(&xTau);
514 if (xTruthParticle)
515 {
516 if (xTruthParticle->isTau())
517 {
518 static const SG::ConstAccessor<char> accIsHadronicTau("IsHadronicTau");
519 if (static_cast<bool>(accIsHadronicTau(*xTruthParticle)))
520 return TruthHadronicTau;
521 else
522 return TruthLeptonicTau;
523 }
524 if (xTruthParticle->isMuon())
525 return TruthMuon;
526 if (xTruthParticle->isElectron())
527 return TruthElectron;
528 }
529
530 // TODO: use const xAOD::Jet* xTruthJet = xAOD::TauHelpers::getLink<xAOD::Jet>(&xTau, "truthJetLink");
531 // currently it is unavailable as templated class is not in icc file
532 static const SG::ConstAccessor< ElementLink< xAOD::JetContainer > > accTruthJetLink("truthJetLink");
533 const ElementLink< xAOD::JetContainer > lTruthParticleLink = accTruthJetLink(xTau);
534 if (lTruthParticleLink.isValid())
535 return TruthJet;
536
537 return Unknown;
538}
539
541{
542 static const SG::ConstAccessor<char> accIsTruthHadronic("IsTruthHadronic");
543 if (!accIsTruthHadronic.isAvailable(xDiTau))
544 Error("TauAnalysisTools::getTruthParticleType", "No truth match information available. Please run DiTauTruthMatchingTool first");
545
546 TruthMatchedParticleType eTruthMatchedParticleType = Unknown;
547
548 if (accIsTruthHadronic(xDiTau))
549 eTruthMatchedParticleType = TruthHadronicDiTau;
550
551 return eTruthMatchedParticleType;
552}
553
554std::vector<const xAOD::TauJet*> TauAnalysisTools::combineTauJetsWithMuonRM(const xAOD::TauJetContainer* taus_std, const xAOD::TauJetContainer* taus_muonRM){
558 // std::string message = "found " + std::to_string(taus_muonRM->size()) + " muon-removal taus";
559 // Info("TauAnalysisTools::getTauJetsWithMuonRM", message.c_str());
560 std::vector<const xAOD::TauJet*> taus_murm_vec(taus_muonRM->begin(), taus_muonRM->end());
561 std::vector<const xAOD::TauJet*> taus_combined;
563 originalTauJetAcc ("originalTauJet");
564 for(const xAOD::TauJet* tau_std : *taus_std){
565 auto replacement_itr = std::find_if(taus_murm_vec.begin(), taus_murm_vec.end(),
566 [&](const xAOD::TauJet* tau_murm){
567 auto link_to_ori_tau = originalTauJetAcc (*tau_murm);
568 if (!link_to_ori_tau.isValid()) { return false; }
569 if (*link_to_ori_tau == tau_std){ return true; }
570 return false;
571 }
572 );
573 if (replacement_itr == taus_murm_vec.end()) { taus_combined.push_back(tau_std); }
574 else {
575 // message = "replacement found at TauJets_MuonRM index " + std::to_string((*replacement_itr)->index()) + " for TauJets index " + std::to_string(tau_std->index());
576 // Info("TauAnalysisTools::getTauJetsWithMuonRM", message.c_str());
577 taus_combined.push_back(*replacement_itr);
578 taus_murm_vec.erase(replacement_itr);
579 // message = std::to_string(taus_murm_vec.size()) + " muon-removal taus left";
580 // Info("TauAnalysisTools::getTauJetsWithMuonRM", message.c_str());
581 }
582 }
583 // Every muon-removal tau should have been used, otherwise there is a problem.
584 assert(taus_murm_vec.empty());
585 return taus_combined;
586}
Scalar deltaR(const MatrixBase< Derived > &vec) const
ATLAS-specific HepMC functions.
static Double_t P(Double_t *tt, Double_t *par)
if(febId1==febId2)
A number of constexpr particle constants to avoid hardcoding them directly in various places.
@ Phi
Definition RPCdef.h:8
@ Eta
Definition RPCdef.h:8
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
Helper class to provide constant type-safe access to aux data.
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition PFO_v1.cxx:95
double etaFinalCalib() const
const PFO * pi0PFO(size_t i) const
Get the pointer to a given pi0 PFO associated with this tau.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition TauJet_v3.cxx:96
virtual double pt() const
The transverse momentum ( ) of the particle.
size_t nPi0PFOs() const
Get the number of pi0 PFO particles associated with this tau.
const TauTrack * track(size_t i, TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged, int *container_index=0) const
Get the pointer to a given tauTrack associated with this tau /*container index needed by trackNonCons...
double ptFinalCalib() const
bool panTauDetail(TauJetParameters::PanTauDetails panTauDetail, int &value) const
Get and set values of pantau details variables via enum.
virtual double eta() const
The pseudorapidity ( ) of the particle.
size_t nTracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
std::vector< const TauTrack * > tracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
Get the v<const pointer> to a given tauTrack collection associated with this tau.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
bool isElectron() const
Whether the particle is an electron (or positron)
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
bool isTau() const
Whether the particle is a tau (or antitau)
bool isMuon() const
Whether the particle is a muon (or antimuon)
static const int PI0
static const int KPLUS
static const int PIPLUS
constexpr double piZeroMassInMeV
the mass of the pi zero (in MeV)
double finalTauEta(const xAOD::TauJet &xTau)
return MVA based tau eta
TruthMatchedParticleType getTruthParticleType(const xAOD::TauJet &xTau)
return TauJet match type
void split(const std::string &sInput, const char cDelim, std::vector< std::string > &vOut)
double finalTauPt(const xAOD::TauJet &xTau)
return MVA based tau pt in GeV
std::vector< const xAOD::TauJet * > combineTauJetsWithMuonRM(const xAOD::TauJetContainer *taus_std, const xAOD::TauJetContainer *taus_muonRM)
combine the standard taujets container with the muon removal container
void correctedPi0Vectors(const xAOD::TauJet *xTau, std::vector< TLorentzVector > &correctedPi0s, TLorentzVector &TauP4)
double finalTauAbsEta(const xAOD::TauJet &xTau)
return MVA based absolute tau eta
void createPi0Vectors(const xAOD::TauJet *xTau, std::vector< TLorentzVector > &vPi0s)
double tauP(const xAOD::TauJet &xTau)
return tau P in GeV
int getNTauDecayParticles(const xAOD::TruthParticle &xTruthTau, int iPdgId, bool bCompareAbsoluteValues)
Count truth matched decay particles of a particular PDGID.
double truthVisTauPt(const xAOD::TauJet &xTau)
return truth match visible tau pt in GeV (if hadronic truth tau match)
double tauPt(const xAOD::TauJet &xTau)
return tau pt in GeV
double tauEta(const xAOD::TauJet &xTau)
return tau eta
double tauLeadTrackEta(const xAOD::TauJet &xTau)
return leading charge tau track eta
const xAOD::TruthParticle * getTruth(const xAOD::TauJet &xTau)
double tauAbsEta(const xAOD::TauJet &xTau)
return absolute tau eta
double truthTauAbsEta(const xAOD::TauJet &xTau)
return truth match tau eta (if hadronic truth tau match)
double finalTauP(const xAOD::TauJet &xTau)
return MVA based tau P in GeV
bool testFileForEOFContainsCharacters(const std::string &sFileName)
returns true if last line in file is empty or the line starts with the number sign
double truthTauPt(const xAOD::TauJet &xTau)
return truth match tau pt in GeV (if hadronic truth tau match)
double truthDecayMode(const xAOD::TauJet &xTau)
return truth decay mode (if hadronic truth tau match)
xAOD::TauJetParameters::DecayMode getTruthDecayMode(const xAOD::TruthParticle &xTruthTau)
Get the Truth Decay Mode from TruthTau particle.
const xAOD::TruthParticle * getTruthParticle(const xAOD::IParticle *, bool debug=false)
return the truthParticle associated to the given IParticle (if any)
PFO_v1 PFO
Definition of the current "pfo version".
Definition PFO.h:17
TauTrack_v1 TauTrack
Definition of the current version.
Definition TauTrack.h:16
TauJet_v3 TauJet
Definition of the current "tau version".
TruthParticle_v1 TruthParticle
Typedef to implementation.
DiTauJet_v1 DiTauJet
Definition of the current version.
Definition DiTauJet.h:17
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".