ATLAS Offline Software
Loading...
Searching...
No Matches
TauCommonCalcVars.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//-----------------------------------------------------------------------------
6// file: TauCommonCalcVars.cxx
7// package: Reconstruction/tauRec
8// authors: Stan Lai
9// date: 2008-05-18
10//
11// This class calculates tau variables after core seed reconstruction
12//-----------------------------------------------------------------------------
14#include <vector>
15
16//-----------------------------------------------------------------------------
17// Constructor
18//-----------------------------------------------------------------------------
19
20TauCommonCalcVars::TauCommonCalcVars(const std::string &name) :
21 TauRecToolBase(name) {
22}
23
24//-----------------------------------------------------------------------------
25// Destructor
26//-----------------------------------------------------------------------------
27
30
31//-----------------------------------------------------------------------------
32// Execution
33//-----------------------------------------------------------------------------
35
37 // Calculate variables that are always valid
39
40 //init some vars
42
43 // Leading track pT and et/pt(lead track)
44 if (pTau.nTracks() > 0) {
45 pTau.setDetail( xAOD::TauJetParameters::leadTrkPt, static_cast<float>( pTau.track(0)->pt() ) );
46
47 float emscale_ptEM = 0.;
48 float emscale_ptHad = 0.;
49
50 if ( !pTau.detail( xAOD::TauJetParameters::etEMAtEMScale, emscale_ptEM ) )
51 {
52 ATH_MSG_DEBUG("retrieval of tau detail failed. stopping calculation of further variables");
53 return StatusCode::SUCCESS;
54 }
55
56 if ( !pTau.detail( xAOD::TauJetParameters::etHadAtEMScale, emscale_ptHad ) )
57 {
58 ATH_MSG_DEBUG("retrieval of tau detail failed. stopping calculation of further variables");
59 return StatusCode::SUCCESS;
60 }
61
62 pTau.setDetail( xAOD::TauJetParameters::etOverPtLeadTrk, static_cast<float>( (emscale_ptEM + emscale_ptHad) / pTau.track(0)->pt() ) );
63 }
64
65 std::vector<const xAOD::TauTrack*> tauTracks = pTau.tracks(xAOD::TauJetParameters::TauTrackFlag::classifiedCharged);
66 for( const xAOD::TauTrack* trk : pTau.tracks((xAOD::TauJetParameters::TauTrackFlag) m_isolationTrackType.value()) ) tauTracks.push_back(trk);
67 if (!tauTracks.empty()) {
68
69 TLorentzVector sumOfTrackVector;
70 double ptSum = 0;
71 double innerPtSum = 0;
72 double sumWeightedDR_tautrack = 0;
73 double innerSumWeightedDR = 0;
74
75 for (const xAOD::TauTrack* tauTrk : tauTracks){
76 sumOfTrackVector += tauTrk->p4();
77
78 double deltaR_tautrack = inTrigger() ? pTau.p4().DeltaR(tauTrk->p4()) : pTau.p4(xAOD::TauJetParameters::IntermediateAxis).DeltaR(tauTrk->p4());
79
80 ptSum += tauTrk->pt();
81 sumWeightedDR_tautrack += deltaR_tautrack * tauTrk->pt();
82
83 //add calculation of innerTrkAvgDist
85 innerPtSum += tauTrk->pt();
86 innerSumWeightedDR += deltaR_tautrack * tauTrk->pt();
87 }
88 }
89 // invariant mass of track system
90 pTau.setDetail( xAOD::TauJetParameters::massTrkSys, static_cast<float>( sumOfTrackVector.M() ) );
91
92 if (ptSum > 0.) {
93 // seedCalo_trkAvgDist
94 pTau.setDetail( xAOD::TauJetParameters::trkAvgDist, static_cast<float>( sumWeightedDR_tautrack / ptSum ) );
95
96 // SumPtTrkFrac
97 pTau.setDetail( xAOD::TauJetParameters::SumPtTrkFrac, static_cast<float>( 1. - innerPtSum/ptSum ) );
98 }
99 else {
102 }
103
104 if (innerPtSum > 0.) {
105 // InnerTrkAvgDist
106 pTau.setDetail( xAOD::TauJetParameters::innerTrkAvgDist, static_cast<float>( innerSumWeightedDR / innerPtSum ) );
107 }
108 else {
110 }
111 }
112
113 return StatusCode::SUCCESS;
114}
#define ATH_MSG_DEBUG(x)
Gaudi::Property< int > m_isolationTrackType
TauCommonCalcVars(const std::string &name="TauCommonCalcVars")
virtual StatusCode execute(xAOD::TauJet &pTau) const override
Execute - called for each tau candidate.
TauRecToolBase(const std::string &name)
bool inTrigger() const
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition TauJet_v3.cxx:96
bool detail(TauJetParameters::Detail detail, int &value) const
Get and set values of common details variables via enum.
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...
void setDetail(TauJetParameters::Detail detail, int value)
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 pt() const
The transverse momentum ( ) of the particle.
TauTrackFlag
Enum for tau track flags.
Definition TauDefs.h:400
@ etHadAtEMScale
Get Hadronic energy at EM scale.
Definition TauDefs.h:196
@ trkAvgDist
Get calibrated EM transverse energy (DEPRECATED since r19)
Definition TauDefs.h:214
@ etEMAtEMScale
Get EM energy at EM scale.
Definition TauDefs.h:194
TauTrack_v1 TauTrack
Definition of the current version.
Definition TauTrack.h:16
TauJet_v3 TauJet
Definition of the current "tau version".