ATLAS Offline Software
Loading...
Searching...
No Matches
TauTruthTrackMatchingTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7// EDM include(s)
11
12using namespace TauAnalysisTools;
13
14//=================================PUBLIC-PART==================================
15//______________________________________________________________________________
16TauTruthTrackMatchingTool::TauTruthTrackMatchingTool( const std::string& name )
17 : AsgTool(name)
18{}
19
20//______________________________________________________________________________
23
24//______________________________________________________________________________
26{
27 return StatusCode::SUCCESS;
28}
29
30//______________________________________________________________________________
31StatusCode TauTruthTrackMatchingTool::classifyTrack(const xAOD::TauTrack& xTrackParticle) const
32{
33 // don't classify tracks if this was already done
34 if (!m_bIsHadronicTrackAvailable.isValid())
35 {
36 static const SG::ConstAccessor<char> accIsHadronicTrack("IsHadronicTrack");
37 bool avail = accIsHadronicTrack.isAvailable(xTrackParticle);
39 if (avail)
40 {
41 ATH_MSG_DEBUG("IsHadronicTrack decoration is available on first track processed, switched of rerun for further taus.");
42 ATH_MSG_DEBUG("If a truth track matching needs to be redone, please pass a shallow copy of the original track.");
43 }
44 }
46 return StatusCode::SUCCESS;
47
48 ATH_CHECK(checkTrackIsTauInheritant(xTrackParticle));
49 ATH_CHECK(checkTrackType(xTrackParticle));
50 return StatusCode::SUCCESS;
51}
52
53//______________________________________________________________________________
54StatusCode TauTruthTrackMatchingTool::classifyTracks(std::vector<const xAOD::TauTrack*>& vTrackParticles) const
55{
56 for (auto xTrackParticle : vTrackParticles)
57 {
58 ATH_CHECK(checkTrackIsTauInheritant(*xTrackParticle));
59 ATH_CHECK(checkTrackType(*xTrackParticle));
60 }
61 return StatusCode::SUCCESS;
62}
63
64//=================================PRIVATE-PART=================================
65//______________________________________________________________________________
66StatusCode TauTruthTrackMatchingTool::checkTrackType(const xAOD::TauTrack& xTrackParticle) const
67{
68 const xAOD::TruthParticle* xTruthParticle = getTruthParticle(xTrackParticle);
69
70 static const SG::Decorator<int> decTruthType("TruthType");
71 if (!xTruthParticle)
72 {
73 decTruthType(xTrackParticle) = TauAnalysisTools::UnclassifiedTrack;
74 return StatusCode::SUCCESS;
75 }
76
77 static const SG::ConstAccessor<float> accTruthMatchProbability("truthMatchProbability");
78
79 if (accTruthMatchProbability(*(xTrackParticle.track())) < 0.5)
80 {
81 decTruthType(xTrackParticle) = TauAnalysisTools::FakeTrack;
82 return StatusCode::SUCCESS;
83 }
84
85 static const SG::ConstAccessor< char > accIsHadronicTrack("IsHadronicTrack");
86 static const SG::ConstAccessor< int > accIsHadronicTrackDecayDepth("IsHadronicTrackDecayDepth");
87 if (static_cast<bool>(accIsHadronicTrack(xTrackParticle)) and accIsHadronicTrackDecayDepth(xTrackParticle) == 0)
88 {
89 decTruthType(xTrackParticle) = TauAnalysisTools::TauTrack;
90 return StatusCode::SUCCESS;
91 }
92
93// AV: The old code had TauAnalysisTools::UnderlyingEventTrack, TauAnalysisTools::PileupTrack, TauAnalysisTools::FakeTrack.
94// Based on https://its.cern.ch/jira/browse/ATLTAU-482 one can speculate the code was designed to distinguish between
95// TauAnalysisTools::UnderlyingEventTrack and TauAnalysisTools::PileupTrack.
96// But such distinction has never been implemented, so the new code has just TauAnalysisTools::UnderlyingEventTrack and TauAnalysisTools::FakeTrack
97
98// protection against code break because of missing uniqueID
99 try{
100 HepMC::no_truth_link(xTruthParticle);
101 } catch (...) {
102 ATH_MSG_DEBUG("unique ID not available for this truth track. Will mark as UnclassifiedTrack ");
103 decTruthType(xTrackParticle) = TauAnalysisTools::UnclassifiedTrack;
104 return StatusCode::SUCCESS;
105 }
106
107 if (HepMC::no_truth_link(xTruthParticle)) decTruthType(xTrackParticle) = TauAnalysisTools::FakeTrack;
108 else if (HepMC::is_simulation_particle(xTruthParticle)) ATH_CHECK(classifyConversion(xTrackParticle, *xTruthParticle));
109 else decTruthType(xTrackParticle) = TauAnalysisTools::UnderlyingEventTrack;
110
111 return StatusCode::SUCCESS;
112}
113
114//______________________________________________________________________________
115StatusCode TauTruthTrackMatchingTool::classifyConversion(const xAOD::TauTrack& xTrackParticle, const xAOD::TruthParticle& xTruthParticle) const
116{
117 static const SG::Decorator<int> decTruthType("TruthType");
118 if (!xTruthParticle.isElectron())
119 {
120 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
121 return StatusCode::SUCCESS;
122 }
123 const xAOD::TruthVertex* xProdVertex = xTruthParticle.prodVtx();
124 if ( !xProdVertex )
125 {
126 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
127 return StatusCode::SUCCESS;
128 }
129 for ( size_t iIncomingParticle = 0; iIncomingParticle < xProdVertex->nIncomingParticles(); ++iIncomingParticle )
130 {
131 const xAOD::TruthParticle* xTruthParent = xProdVertex->incomingParticle(iIncomingParticle);
132 if (!xTruthParent)
133 {
134 ATH_MSG_WARNING("Truth parent of tau decay was not found in TruthParticles container. Please ensure that this container has the full tau decay information or produce the TruthTaus container in AtlasDerivation.\nInformation on how to do this can be found here:\nhttps://twiki.cern.ch/twiki/bin/viewauth/AtlasProtected/TauPreRecommendations2015#Accessing_Tau_Truth_Information");
135 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
136 return StatusCode::SUCCESS;
137 }
138
139 if (!xTruthParent->isPhoton())
140 {
141 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
142 return StatusCode::SUCCESS;
143 }
144 }
145
146 size_t iElectrons = 0;
147 for ( size_t iOutgoingParticle = 0; iOutgoingParticle < xProdVertex->nOutgoingParticles(); ++iOutgoingParticle )
148 {
149 const xAOD::TruthParticle* xTruthDaughter = xProdVertex->outgoingParticle(iOutgoingParticle);
150 if (!xTruthDaughter)
151 {
152 ATH_MSG_WARNING("Truth daughter of tau decay was not found in TruthParticles container. Please ensure that this container has the full tau decay information or produce the TruthTaus container in AtlasDerivation.\nInformation on how to do this can be found here:\nhttps://twiki.cern.ch/twiki/bin/viewauth/AtlasProtected/TauPreRecommendations2015#Accessing_Tau_Truth_Information");
153 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
154 return StatusCode::SUCCESS;
155 }
156
157 if (!xTruthDaughter->isElectron())
158 {
159 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
160 return StatusCode::SUCCESS;
161 }
162 iElectrons++;
163 }
164 if (iElectrons != 2)
165 {
166 decTruthType(xTrackParticle) = TauAnalysisTools::SecondaryTrack;
167 return StatusCode::SUCCESS;
168 }
169 decTruthType(xTrackParticle) = TauAnalysisTools::ConversionTrack;
170 return StatusCode::SUCCESS;
171}
172
173//______________________________________________________________________________
175{
176 static const SG::ConstAccessor< ElementLink<xAOD::TruthParticleContainer> > accTruthParticleLink("truthParticleLink");
177 auto xTruthParticleLink = accTruthParticleLink(*(xTrackParticle.track()));
178 //check validity of truth particle element link
179 if (xTruthParticleLink.isValid())
180 return (*xTruthParticleLink);
181 return nullptr;
182}
183
184//______________________________________________________________________________
186{
187 static const SG::Decorator< char > decIsHadronicTrack("IsHadronicTrack");
188 static const SG::Decorator< int > decIsHadronicTrackDecayDepth("IsHadronicTrackDecayDepth");
189 static const SG::Decorator< std::string > decDecayHistory("DecayHistory");
190 decIsHadronicTrack(xTrackParticle) = (char)false;
191 int iDepth = -1;
192 const xAOD::TruthParticle* xTruthParticle = getTruthParticle(xTrackParticle);
193 std::string sHistory = "";
194 if (xTruthParticle)
195 {
196 sHistory = std::to_string(xTruthParticle->pdgId());
197 if (checkTruthParent(*xTruthParticle, iDepth, sHistory))
198 {
199 decIsHadronicTrack(xTrackParticle) = (char)true;
200 }
201 }
202 decIsHadronicTrackDecayDepth(xTrackParticle) = iDepth;
203 decDecayHistory(xTrackParticle) = sHistory;
204
205 return StatusCode::SUCCESS;
206}
207
208//______________________________________________________________________________
209bool TauTruthTrackMatchingTool::checkTruthParent(const xAOD::TruthParticle& xTruthParticle, int& iDepth, std::string& sHistory) const
210{
211 iDepth++;
212 if (xTruthParticle.hasProdVtx())
213 {
214 const xAOD::TruthVertex* xVertex = xTruthParticle.prodVtx();
215 if (xVertex->nIncomingParticles() != 1)
216 ATH_MSG_VERBOSE("not exactly one incomming particles for vertex, number of incomming particles: "<<xVertex->nIncomingParticles());
217 if (xVertex->nIncomingParticles() > 0)
218 {
219 const xAOD::TruthParticle* xTruthParticleParent = xVertex->incomingParticle(0);
220 if (xTruthParticleParent) {
221 // store parent pdgID in history
222 sHistory.insert(0, std::to_string(xTruthParticleParent->pdgId())+":");//xTruthParticleParent->pdgId());
223 if (MC::isTau(xTruthParticleParent))
224 {
225 return true;
226 }
227 else
228 {
229 return checkTruthParent(*xTruthParticleParent, iDepth, sHistory);
230 }
231 } else {
232 ATH_MSG_WARNING("vertex has incoming particles but no valid parent particle");
233 }
234 }
235 }
236 return false;
237}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
ATLAS-specific HepMC functions.
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.
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
const xAOD::TruthParticle * getTruthParticle(const xAOD::TauTrack &xTrackParticle) const
virtual StatusCode classifyTrack(const xAOD::TauTrack &xTrackParticle) const override
CxxUtils::CachedValue< bool > m_bIsHadronicTrackAvailable
virtual StatusCode classifyTracks(std::vector< const xAOD::TauTrack * > &vTracks) const override
StatusCode checkTrackIsTauInheritant(const xAOD::TauTrack &xTrackParticle) const
bool checkTruthParent(const xAOD::TruthParticle &xTruthParticle, int &iDepth, std::string &sHistory) const
StatusCode classifyConversion(const xAOD::TauTrack &xTrackParticle, const xAOD::TruthParticle &xTruthParticle) const
virtual StatusCode initialize() override
Declare the interface that the class provides.
virtual ASG_TOOL_CLASS(TauTruthTrackMatchingTool, TauAnalysisTools::ITauTruthTrackMatchingTool) public ~TauTruthTrackMatchingTool()
Create a proper constructor for Athena.
StatusCode checkTrackType(const xAOD::TauTrack &xTrackParticle) const
const TrackParticle * track() const
int pdgId() const
PDG ID code.
bool hasProdVtx() const
Check for a production vertex on this particle.
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
bool isElectron() const
Whether the particle is an electron (or positron)
bool isPhoton() const
Whether the particle is a photon.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
size_t nIncomingParticles() const
Get the number of incoming particles.
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
bool no_truth_link(const T &p)
Method to establish if a if the object is linked to something which was never saved to the HepMC Trut...
bool isTau(const T &p)
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TauTrack_v1 TauTrack
Definition of the current version.
Definition TauTrack.h:16
TruthParticle_v1 TruthParticle
Typedef to implementation.