ATLAS Offline Software
Loading...
Searching...
No Matches
InDetTrackTruthOriginTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
7
12
13#include <math.h>
14
15namespace InDet {
16
17 InDetTrackTruthOriginTool::InDetTrackTruthOriginTool(const std::string& name) :
18 asg::AsgTool(name)
19 {
20
21#ifndef XAOD_STANDALONE
22 declareInterface<IInDetTrackTruthOriginTool>(this);
23#endif
24
25 declareProperty("matchingProbabilityCut", m_matchingProbabilityCut = 0.5);
26
27 declareProperty("truthParticleLinkName", m_truthParticleLinkName = "truthParticleLink");
28 declareProperty("truthMatchProbabilityAuxName", m_truthMatchProbabilityAuxName = "truthMatchProbability");
29 declareProperty("isFullPileUpTruth", m_isFullPileupTruth = false);
30
31 }
32
33 InDetTrackTruthOriginTool::~InDetTrackTruthOriginTool() = default;
34
36 if(m_isFullPileupTruth) ATH_MSG_INFO("InDetTrackTruthOriginTool configured for full pile-up truth sample");
37 return StatusCode::SUCCESS;
38 }
39
41
42 // if the track doesnt't have a valid truth link, skip to the next track
43 // in practice, all tracks seem to have a truth link, but we need to also
44 // check whether it's valid
47 truthParticleLinkAcc ("truthParticleLink");
48 if ( !truthParticleLinkAcc.isAvailable(*track) ) {
49 return nullptr;
50 }
51
52 // retrieve the link and check its validity
53 const TruthLink &link = truthParticleLinkAcc(*track);
54
55 // a missing or invalid link implies truth particle has been dropped from
56 // the truth record at some stage - probably it was from pilup which by
57 // default we don't store truth information for
58 if(!link or !link.isValid()) {
59 return nullptr;
60 }
61
62 // seems safe to access and return the linked truth particle
63 return *link;
64 }
65
67
68 int origin = 0;
69
70 // is it fragmentation? yes until proven otherwise
71 bool isFragmentation = true;
72
73 // from B decay chain?
74 if(isFrom(truth, MC::BQUARK)) {
75 origin = origin | (0x1 << InDet::TrkOrigin::BHadronDecay);
76 isFragmentation = false;
77 }
78
79 // from D decay chain?
80 if(isFrom(truth, MC::CQUARK)) {
81 origin = origin | (0x1 << InDet::TrkOrigin::DHadronDecay);
82 isFragmentation = false;
83 }
84
85 // from tau decay chain?
86 if(isFrom(truth, MC::TAU)) {
87 origin = origin | (0x1 << InDet::TrkOrigin::TauDecay);
88 isFragmentation = false;
89 }
90
91 // Secondary?
93
94 // unknown origin (e.g. parent not in the record)
95 if(truth->nParents() != 1) {
96 origin = origin | (0x1 << InDet::TrkOrigin::OtherOrigin);
97 }
98
99 else {
100 const xAOD::TruthParticle *parent = truth->parent(0);
101
102 if(parent == nullptr) {
103 origin = origin | (0x1 << InDet::TrkOrigin::OtherOrigin);
104 }
105
106 // sub-categorize secondaries which have one valid parent
107 else {
108
109 int pdgId = truth->pdgId();
110 int parentId = parent->pdgId();
111
112 // photon conversion
113 if(parent->isPhoton() && truth->isElectron()) {
114 origin = origin | (0x1 << InDet::TrkOrigin::GammaConversion);
115 }
116
117 // Strange Mesons
118 else if(parent->isStrangeMeson() && parent->nChildren() == 2) {
119 origin = origin | (0x1 << InDet::TrkOrigin::StrangeMesonDecay);
120 // specifically Kshort
121 if (abs(pdgId) == 211 && parentId == 310) {
122 origin = origin | (0x1 << InDet::TrkOrigin::KshortDecay);
123 }
124 }
125
126 // Strange Baryons
127 else if(parent->isStrangeBaryon() && parent->nChildren() == 2) {
128 origin = origin | (0x1 << InDet::TrkOrigin::StrangeBaryonDecay);
129 // specifically Lambdas
130 if ((abs(pdgId) == MC::PIPLUS || abs(pdgId) == MC::PROTON) && abs(parentId) == MC::LAMBDA0) {
131 origin = origin | (0x1 << InDet::TrkOrigin::LambdaDecay);
132 }
133 }
134
135 // other long living particle decays
136 else if(parent->isHadron() && parent->nChildren() == 2) {
137 origin = origin | (0x1 << InDet::TrkOrigin::OtherDecay);
138 }
139
140 // hadronic interactions
141 else if(parent->nChildren() > 2) {
142 origin = origin | (0x1 << InDet::TrkOrigin::HadronicInteraction);
143 }
144
145 // other secondaries
146 else {
147 origin = origin | (0x1 << InDet::TrkOrigin::OtherSecondary);
148 }
149 }
150 }
151
152 isFragmentation = false;
153 }
154
155 // uncategorized: call it fragmentation
156 if(isFragmentation) {
157 origin = origin | (0x1 << InDet::TrkOrigin::Fragmentation);
158 }
159
160 return origin;
161 }
162
164
165 // get truth particle
166 const xAOD::TruthParticle* truth = getTruth( track );
167
168 // get track TMP
170 float truthProb = tmpAcc( *track );
171
172 int origin = 0;
173
174 // truth link is broken: call it from pileup (not true for < 100 MeV and some specialised samples!)
175 if (!truth){
176 origin = origin | (0x1 << InDet::TrkOrigin::Pileup);
177 }
178
179 // For full pile-up truth samples, check explicitly if truth record is present among hard-scatter collection to assign Pileup origin
180 else if(m_isFullPileupTruth){
181 const xAOD::TruthEventContainer* truthEventContainer(nullptr);
182 if(evtStore()->retrieve(truthEventContainer, "TruthEvents").isFailure()){
183 ATH_MSG_ERROR("InDetTrackTruthOriginTool configured for full pile-up truth but could not retrieve TruthEvents container");
184 }
185 const xAOD::TruthEvent* event = truthEventContainer ? truthEventContainer->at(0) : nullptr;
186
187 if(event){
188 const auto& links = event->truthParticleLinks();
189
190 bool isFromHSProdVtx = false;
191 for (const auto& link : links){
192 if(link.isValid() && truth == *link){
193 isFromHSProdVtx = true;
194 break;
195 }
196 }
197
198 if(!isFromHSProdVtx) origin = origin | (0x1 << InDet::TrkOrigin::Pileup);
199 }
200 }
201
202 // low TruthMatchProbability: call it fake (also includes poorly reconstructed tracks)
203 if(truthProb < m_matchingProbabilityCut) {
204 origin = origin | (0x1 << InDet::TrkOrigin::Fake);
205 }
206
207 // truth link is present: find truth origin
208 if (truth) {
209 origin = origin | getTruthOrigin(truth);
210 }
211
212 return origin;
213 }
214
215 bool InDetTrackTruthOriginTool::isFrom(const xAOD::TruthParticle* truth, int flav) const {
216 return isFromRec( truth, flav, 0 );
217 }
218
219 bool InDetTrackTruthOriginTool::isFromRec(const xAOD::TruthParticle* truth, int flav, int depth) const {
220
221 if ( truth == nullptr ) return false;
222
223 if ( depth > 30 ) return false;
224
225 if( flav != MC::BQUARK && flav != MC::CQUARK && flav != MC::TAU ) return false;
226
227 if( flav == MC::BQUARK && truth->isBottomHadron() ) return true;
228
229 if( flav == MC::CQUARK && truth->isCharmHadron() ) return true;
230
231 if( flav == MC::TAU && MC::isTau(truth) ) return true;
232
233
234 for(unsigned int p=0; p<truth->nParents(); p++) {
235 const xAOD::TruthParticle* parent = truth->parent(p);
236 if(parent == truth ) continue ; // avoid infinite recursion
237 if( isFromRec(parent, flav, depth+1) ) return true;
238 }
239
240 return false;
241 }
242
243}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
Helper class to provide constant type-safe access to aux data.
ATLAS-specific HepMC functions.
ElementLink< xAOD::TruthParticleContainer > TruthLink
ServiceHandle< StoreGateSvc > & evtStore()
const T * at(size_type n) const
Access an element, as an rvalue.
virtual bool isFromRec(const xAOD::TruthParticle *truth, int flav, int depth=0) const
recursion-safe(r) version of isFrom
virtual int getTruthOrigin(const xAOD::TruthParticle *truth) const override
Computes the truth particle origin.
virtual StatusCode initialize() override
virtual const xAOD::TruthParticle * getTruth(const xAOD::TrackParticle *track) const override
Safely access a track's linked truth particle, if available.
virtual int getTrackOrigin(const xAOD::TrackParticle *track) const override
Computes the track origin.
virtual bool isFrom(const xAOD::TruthParticle *truth, int flav) const override
Check if a truth particle is from the specified origin (from B or D hadron, or tau)
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:569
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.
bool isBottomHadron() const
Determine if the PID is that of a b-hadron.
int pdgId() const
PDG ID code.
const TruthParticle_v1 * parent(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
size_t nParents() const
Number of parents of this particle.
bool isElectron() const
Whether the particle is an electron (or positron)
bool isCharmHadron() const
Determine if the PID is that of a c-hadron.
std::string depth
tag string for intendation
Definition fastadd.cxx:46
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...
Primary Vertex Finder.
static const int CQUARK
static const int TAU
static const int PIPLUS
static const int BQUARK
static const int LAMBDA0
bool isTau(const T &p)
static const int PROTON
TruthEventContainer_v1 TruthEventContainer
Declare the latest version of the truth event container.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthEvent_v1 TruthEvent
Typedef to implementation.
Definition TruthEvent.h:17
TruthParticle_v1 TruthParticle
Typedef to implementation.