ATLAS Offline Software
Loading...
Searching...
No Matches
RecoToTruthAssociationAlg.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
11#include "TrkTrack/Track.h"
18using namespace xAOD::P4Helpers;
19namespace {
20 constexpr unsigned int dummy_unsigned = 999;
21 constexpr int com_bit = (1<<xAOD::Muon::Author::Commissioning);
22 void increment_unsigned(unsigned& val) {
23 if (val == dummy_unsigned)
24 val = 1;
25 else
26 ++val;
27 }
29 const SG::ConstAccessor<int> acc_origin("truthOrigin");
30 const SG::ConstAccessor<int> acc_type("truthType");
31 const SG::ConstAccessor<TruthLink_t> acc_link("truthParticleLink");
32 //
33 const SG::Decorator<int> dec_origin("truthOrigin");
34 const SG::Decorator<int> dec_type ("truthType");
35 const SG::Decorator<TruthLink_t> dec_link("truthParticleLink");
36
37 static const SG::ConstAccessor<std::vector<unsigned long long>> truthMdtHitsAcc("truthMdtHits");
38 static const SG::ConstAccessor<std::vector<unsigned long long>> truthCscHitsAcc("truthCscHits");
39 static const SG::ConstAccessor<std::vector<unsigned long long>> truthRpcHitsAcc("truthRpcHits");
40 static const SG::ConstAccessor<std::vector<unsigned long long>> truthTgcHitsAcc("truthTgcHits");
41
42} // namespace
43
44namespace Muon {
45using namespace MuonStationIndex;
46// Initialize method:
48 ATH_CHECK(m_idHelperSvc.retrieve());
49 ATH_CHECK(m_truthMuKey.initialize());
50 ATH_CHECK(m_recoMuKey.initialize());
51
52 if (m_recoLink.empty()){
54 } else {
56 }
57
58 ATH_CHECK(m_muonTruthRecoLink.initialize(!m_recoLink.empty()));
65 for (const std::string& trk_coll : m_assocTrkContainers.value()){
66 m_inputDecorKey.emplace_back(trk_coll + ".truthParticleLink");
67 }
68
69 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(acc_origin.auxid()));
70 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(acc_type.auxid()));
71 if (m_idHelperSvc->hasMDT()) {
72 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(truthMdtHitsAcc.auxid()));
73 }
74 if (m_idHelperSvc->hasRPC()) {
75 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(truthRpcHitsAcc.auxid()));
76 }
77 if (m_idHelperSvc->hasTGC()) {
78 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(truthTgcHitsAcc.auxid()));
79 }
80 if (m_idHelperSvc->hasCSC()) {
81 m_inputDecorKey.emplace_back(m_truthMuKey, SG::AuxTypeRegistry::instance().getName(truthCscHitsAcc.auxid()));
82 }
83
84 ATH_CHECK(m_inputDecorKey.initialize());
85 return StatusCode::SUCCESS;
86}
87
88// Execute method:
89StatusCode RecoToTruthAssociationAlg::execute(const EventContext& ctx) const {
90 const xAOD::TruthParticleContainer* muonTruthContainer{nullptr};
91 ATH_CHECK(SG::get(muonTruthContainer, m_truthMuKey , ctx));
92 std::unique_ptr<SG::WriteDecorHandle<xAOD::TruthParticleContainer, ElementLink<xAOD::MuonContainer>>> muonTruthParticleRecoLink{};
93 if (!m_muonTruthRecoLink.empty()) {
94 muonTruthParticleRecoLink = std::make_unique<SG::WriteDecorHandle<xAOD::TruthParticleContainer, ElementLink<xAOD::MuonContainer>>>(m_muonTruthRecoLink, ctx);
95 }
96
98 ctx);
99 if (!muonTruthParticleLink.isValid()) {
100 ATH_MSG_WARNING("muon particle container not valid");
101 return StatusCode::FAILURE;
102 }
106 ctx);
109
110 // add link to reco muons and viceversa
111 bool saw_staco = false;
112 bool decor_staco = false;
113
114 // loop over muons
115 for (const xAOD::Muon* muon : *muonTruthParticleLink) {
116 // use primary track particle to get the truth link (except for the case of STACO, where we must use the ID track particle, as the
117 // combined is not truth-matched)
118 ATH_MSG_DEBUG("muon with pT " << muon->pt() << " MeV, eta: " << muon->eta() << ", phi " << muon->phi() << " and author "
119 << muon->author());
120 const xAOD::TrackParticle* tp = nullptr;
121 if (m_associateWithInDetTP || muon->author() == xAOD::Muon::STACO || muon->author() == xAOD::Muon::MuGirl) {
122 tp = muon->trackParticle(xAOD::Muon::InnerDetectorTrackParticle);
123 } else {
124 tp = muon->primaryTrackParticle();
125 }
126
127 bool foundTruth{false}, setOrigin{false};
128 if (tp) {
129 // Associate reco with truth muon. Loop over reconstructed muons, get track particle for each one.
130 // Each track particle should carry a link to the corresponding truth particle. Then compare this truth particle link with the
131 // given truth muon particle
132 if (acc_origin.isAvailable(*tp) && acc_origin(*tp) != 0) {
133 muonTruthParticleOrigin(*muon) = acc_origin(*tp);
134 muonTruthParticleType(*muon) = acc_type(*tp);
135 setOrigin = true;
136 }
137
139 if (acc_link.isAvailable(*tp)) {
140 truthLink = acc_link(*tp);
141 } else {
142 ATH_MSG_DEBUG("Could not find any truth link associated with track having pt:"<<tp->pt()<<" MeV, eta: "<<tp->eta()<<", phi: "<<tp->phi()<<", charge: "<<tp->charge()<<". d0:"<<tp->d0()<<", z0: "<<tp->z0());
143 }
144
145 if (truthLink.isValid()) {
146 ATH_MSG_VERBOSE(" Got valid truth link for muon author " << muon->author() << " uniqueID " << HepMC::uniqueID(*truthLink));
147 // loop over truth particles
148 for (const xAOD::TruthParticle* truthParticle : *muonTruthContainer) {
149 if (!MC::isStable(truthParticle)) continue;
150 ATH_MSG_DEBUG("Got truth muon with uniqueID " << HepMC::uniqueID(truthParticle) << " pt " << truthParticle->pt());
151 if ( !HepMC::is_sim_descendant(*truthLink, truthParticle)) {
152 ATH_MSG_VERBOSE("UniqueID truth link: " << HepMC::uniqueID(*truthLink)
153 << " is not decendant of " << HepMC::uniqueID(truthParticle));
154 continue;
155 }
156 ATH_MSG_VERBOSE("Truth muon uniqueID matches -> creating link with truth particle " << HepMC::uniqueID(*truthLink));
157 foundTruth = true;
159 ElementLink<xAOD::TruthParticleContainer> muonTruthLink{*muonTruthContainer,
160 truthParticle->index(),
161 ctx};
162 muonTruthLink.toPersistent();
163 muonTruthParticleLink(*muon) = muonTruthLink;
164 if (!setOrigin) {
165 muonTruthParticleOrigin(*muon) = acc_origin(*tp);
166 muonTruthParticleType(*muon) = acc_type(*tp);
167 setOrigin = true;
168 }
170 if (muonTruthParticleRecoLink && muonTruthParticleRecoLink->operator()(*truthParticle).isValid()){
171 const xAOD::Muon* decor_muon = *muonTruthParticleRecoLink->operator()(*truthParticle);
172 ATH_MSG_VERBOSE("Truth particle is already decorated with reco muon "<<decor_muon->pt()*1.e-3
173 <<" eta: "<<decor_muon->eta()<<" phi: "<<decor_muon->phi()<<" charge: "<<
174 decor_muon->charge()<<" author: "<<decor_muon->author()<<" all authors: "<<
175 decor_muon->allAuthors());
176
177 // Check first if the exiting muon has a better author
178 if (MuonCombined::authorRank(decor_muon->author()) < MuonCombined::authorRank(muon->author())){
179 ATH_MSG_DEBUG("Author of the decorated muon is better than the one of the new candidate");
180 continue;
181 }
183 const int com_score = (muon->allAuthors() & com_bit) - (decor_muon->allAuthors() &com_bit);
184 if (com_score > 0){
185 ATH_MSG_DEBUG("Found two muons reconstructed by an equivalent author. But this one is from the commissioning chain");
186 continue;
187 }
189 if (deltaR2(muon,truthParticle) >= deltaR2(muon, decor_muon)) continue;
190 }
191
192
193 ElementLink<xAOD::MuonContainer> muonLink{muon, *muonTruthParticleLink, ctx};
194
196 std::vector<unsigned int> nprecHitsPerChamberLayer(toInt(ChIndex::ChIndexMax), dummy_unsigned);
197 std::vector<unsigned int> nphiHitsPerChamberLayer(toInt(PhiIndex::PhiIndexMax), dummy_unsigned);
198 std::vector<unsigned int> ntrigEtaHitsPerChamberLayer(toInt(PhiIndex::PhiIndexMax), dummy_unsigned);
199
200 constexpr int author_sel = (1<<xAOD::Muon::MuidCo) | (1<<xAOD::Muon::MuidSA) | (1<<xAOD::Muon::MuGirl);
201 count_chamber_layers(muon->allAuthors() & author_sel
202 ? truthParticle
203 : nullptr,
204 tp->track(), nprecHitsPerChamberLayer, nphiHitsPerChamberLayer, ntrigEtaHitsPerChamberLayer);
206 muonTruthParticleNPrecMatched(*muon) = nprecHitsPerChamberLayer;
207 muonTruthParticleNPhiMatched(*muon) = nphiHitsPerChamberLayer;
208 muonTruthParticleNTrigEtaMatched(*muon) = ntrigEtaHitsPerChamberLayer;
209
210 if (muonTruthParticleRecoLink) (*muonTruthParticleRecoLink)(*truthParticle) = muonLink;
211 break;
212 }
213 } else {
214 ATH_MSG_DEBUG("Invalid truth link");
215 }
216 } else {
217 ATH_MSG_WARNING("Could not find the appropiate track particle for muon with pT: " << muon->pt() * 1.e-3 << " GeV, eta: "
218 << muon->eta() << ", phi: " << muon->phi()
219 << " author: " << muon->author());
220 }
221
222 if (!setOrigin) {
223 muonTruthParticleOrigin(*muon) = 0;
224 muonTruthParticleType(*muon) = 0;
225 }
226 if (!foundTruth) {
227 muonTruthParticleLink(*muon) = ElementLink<xAOD::TruthParticleContainer>();
228 // add these empty vectors
229 muonTruthParticleNPrecMatched(*muon) = std::vector<unsigned int>{};
230 muonTruthParticleNPhiMatched(*muon) = std::vector<unsigned int>{};
231 muonTruthParticleNTrigEtaMatched(*muon) = std::vector<unsigned int>{};
232 }
235 if (muon->author() == xAOD::Muon::STACO) {
236 const xAOD::TrackParticle* cmb_trk = muon->trackParticle(xAOD::Muon::CombinedTrackParticle);
237 if (!cmb_trk){
238 ATH_MSG_WARNING("Even a STACO muon should have a combined track");
239 continue;
240 } else {
241 if (!saw_staco) {
242 saw_staco = true;
243 decor_staco = !dec_origin.isAvailable (*cmb_trk);
244 }
245 if (decor_staco) {
246 dec_origin(*cmb_trk) = acc_origin(*muon);
247 dec_type(*cmb_trk) = acc_type(*muon);
248 dec_link(*cmb_trk) = acc_link(*muon);
249 }
250 }
251 }
252 }
254 if (muonTruthParticleRecoLink && !muonTruthParticleRecoLink->isAvailable()) {
255 for (const xAOD::TruthParticle* truthParticle : **muonTruthParticleRecoLink) {
256 ATH_MSG_DEBUG("no reco muon link set, add an empty one");
257 (*muonTruthParticleRecoLink)(*truthParticle) = ElementLink<xAOD::MuonContainer>();
258 }
259 }
260
261 return StatusCode::SUCCESS;
262}
263
265 std::vector<unsigned int>& nprecHitsPerChamberLayer,
266 std::vector<unsigned int>& nphiHitsPerChamberLayer,
267 std::vector<unsigned int>& ntrigEtaHitsPerChamberLayer) const {
268
269 if (!truthParticle || !truthMdtHitsAcc.isAvailable(*truthParticle)) {
270 ATH_MSG_DEBUG("muon has no truth hits vector in the truth association alg");
271 nprecHitsPerChamberLayer.clear();
272 nphiHitsPerChamberLayer.clear();
273 ntrigEtaHitsPerChamberLayer.clear();
274 return;
275 }
276 const std::vector<unsigned long long>& mdtTruth = truthMdtHitsAcc(*truthParticle);
277 std::vector<unsigned long long> cscTruth;
278
279 if (m_idHelperSvc->hasCSC()) cscTruth = truthCscHitsAcc(*truthParticle);
280 const std::vector<unsigned long long>& rpcTruth = truthRpcHitsAcc(*truthParticle);
281 const std::vector<unsigned long long>& tgcTruth = truthTgcHitsAcc(*truthParticle);
282
283 for (const Trk::TrackStateOnSurface* tsit : *ptrk->trackStateOnSurfaces()) {
284 if (!tsit || !tsit->trackParameters() || !tsit->measurementOnTrack()) continue;
285 const Trk::MeasurementBase* meas = tsit->measurementOnTrack();
286 Identifier id;
287 const Trk::RIO_OnTrack* rot = dynamic_cast<const Trk::RIO_OnTrack*>(meas);
288 if (rot)
289 id = rot->identify();
290 else {
291 const Muon::CompetingMuonClustersOnTrack* crot = dynamic_cast<const Muon::CompetingMuonClustersOnTrack*>(meas);
292 if (crot && !crot->containedROTs().empty() && crot->containedROTs().front()) id = crot->containedROTs().front()->identify();
293 }
294 if (!m_idHelperSvc->isMuon(id)) continue;
295
296 bool measPhi = m_idHelperSvc->measuresPhi(id);
297 bool isTgc = m_idHelperSvc->isTgc(id);
298 ChIndex chIndex = !isTgc ? m_idHelperSvc->chamberIndex(id) : ChIndex::ChUnknown;
299 if (m_idHelperSvc->isMdt(id)) {
300 for (unsigned int i = 0; i < mdtTruth.size(); ++i) {
301 if (id == mdtTruth[i]) {
302 if (chIndex != ChIndex::ChUnknown) {
303 increment_unsigned(nprecHitsPerChamberLayer.at(toInt(chIndex)));
304 }
305 break;
306 }
307 }
308 } else if (m_idHelperSvc->hasCSC() && m_idHelperSvc->isCsc(id)) {
309 for (unsigned int i = 0; i < cscTruth.size(); ++i) {
310 if (id != cscTruth[i]) continue;
311 if (measPhi) {
312 PhiIndex index = m_idHelperSvc->phiIndex(id);
313 if (index != PhiIndex::PhiUnknown) {
314 increment_unsigned(nphiHitsPerChamberLayer.at(toInt(index)));
315 }
316 } else {
317 if (chIndex != ChIndex::ChUnknown) {
318 increment_unsigned(nprecHitsPerChamberLayer.at(toInt(chIndex)));
319 }
320 }
321 break;
322 }
323 } else if (m_idHelperSvc->isRpc(id)) {
324 for (unsigned int i = 0; i < rpcTruth.size(); ++i) {
325 if (id != rpcTruth[i]) { continue; }
326 PhiIndex index = m_idHelperSvc->phiIndex(id);
327 if (index != PhiIndex::PhiUnknown) {
328 if (measPhi) {
329 increment_unsigned(nphiHitsPerChamberLayer.at(toInt(index)));
330 } else {
331 increment_unsigned(ntrigEtaHitsPerChamberLayer.at(toInt(index)));
332 }
333 }
334 break;
335 }
336 } else if (m_idHelperSvc->isTgc(id)) {
337 for (unsigned int i = 0; i < tgcTruth.size(); ++i) {
338 if (id != tgcTruth[i]) { continue; }
339 PhiIndex index = m_idHelperSvc->phiIndex(id);
340 if (index != PhiIndex::PhiUnknown) {
341 if (measPhi) {
342 increment_unsigned(nphiHitsPerChamberLayer.at(toInt(index)));
343 } else {
344 increment_unsigned(ntrigEtaHitsPerChamberLayer.at(toInt(index)));
345 }
346 }
347 break;
348 }
349 }
350 } // end loop over TSOS
351 ATH_MSG_DEBUG("finished loop over TSOS");
352
353 // now, have to check if there are non-zero truth hits in indices without reco hits
354 clear_dummys(mdtTruth, nprecHitsPerChamberLayer);
355 clear_dummys(cscTruth, nprecHitsPerChamberLayer);
356
357 clear_dummys(cscTruth, nphiHitsPerChamberLayer);
358 clear_dummys(rpcTruth, nphiHitsPerChamberLayer);
359 clear_dummys(tgcTruth, nphiHitsPerChamberLayer);
360
361 clear_dummys(rpcTruth, ntrigEtaHitsPerChamberLayer);
362 clear_dummys(tgcTruth, ntrigEtaHitsPerChamberLayer);
363}
364void RecoToTruthAssociationAlg::clear_dummys(const std::vector<unsigned long long>& identifiers, std::vector<unsigned int>& vec) const {
367 if (identifiers.empty()) { return; }
368 for (unsigned int i = 0; i < vec.size(); ++i) {
369 if (vec[i] != dummy_unsigned) continue;
370 for (unsigned j = 0; j < identifiers.size(); ++j) {
371 const Identifier id{identifiers[j]};
372 if (m_idHelperSvc->measuresPhi(id)) {
373 const auto phiIdx = static_cast<PhiIndex>(i);
374 if (m_idHelperSvc->phiIndex(id) == phiIdx) {
375 vec[i] = 0;
376 break;
377 }
378 } else {
379 const auto chIdx = static_cast<ChIndex>(i);
380 if (m_idHelperSvc->chamberIndex(id) == chIdx) {
381 vec[i] = 0;
382 break;
383 }
384 }
385 }
386 }
387}
388}
#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)
std::vector< size_t > vec
Helper class to provide constant type-safe access to aux data.
ATLAS-specific HepMC functions.
Handle class for adding a decoration to an object.
ElementLink< xAOD::TruthParticleContainer > TruthLink_t
Class for competing MuonClusters, it extends the Trk::CompetingRIOsOnTrack base class.
const std::vector< const MuonClusterOnTrack * > & containedROTs() const
returns the vector of SCT_ClusterOnTrack objects .
SG::ReadHandleKey< xAOD::MuonContainer > m_recoMuKey
SG::ReadDecorHandleKeyArray< xAOD::IParticleContainer > m_inputDecorKey
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleOrigin
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleNPhiMatched
Gaudi::Property< bool > m_associateWithInDetTP
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_muonTruthRecoLink
Gaudi::Property< std::vector< std::string > > m_assocTrkContainers
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleType
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleLink
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleNTrigEtaMatched
void count_chamber_layers(const xAOD::IParticle *truth_particle, const Trk::Track *ptrk, std::vector< unsigned int > &nprecHitsPerChamberLayer, std::vector< unsigned int > &nphiHitsPerChamberLayer, std::vector< unsigned int > &ntrigEtaHitsPerChamberLayer) const
void clear_dummys(const std::vector< unsigned long long > &identifiers, std::vector< unsigned int > &vec) const
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
Gaudi::Property< std::string > m_recoLink
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthMuKey
Key to the filtered muon truth particles.
StatusCode execute(const EventContext &ctx) const override
SG::WriteDecorHandleKey< xAOD::MuonContainer > m_muonTruthParticleNPrecMatched
static AuxTypeRegistry & instance()
Return the singleton registry instance.
Helper class to provide constant type-safe access to aux data.
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
Handle class for adding a decoration to an object.
This class is the pure abstract base class for all fittable tracking measurements.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Identifier identify() const
return the identifier -extends MeasurementBase
represents the track state (measurement, material, fit parameters and quality) at a surface.
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
Class providing the definition of the 4-vector interface.
XAOD_AUXDATA_DEPRECATED bool isAvailable(const std::string &name, const std::string &clsname="") const
Check if a user property is available for reading or not.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
float charge() const
uint16_t allAuthors() const
Get all the authors of this Muon.
Author author() const
int uniqueID(const T &p)
bool is_sim_descendant(const T1 &p1, const T2 &p2)
Method to check if the first particle is a descendant of the second in the simulation,...
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
int authorRank(const xAOD::Muon::Author &a)
Definition TagBase.h:23
ChIndex chIndex(const std::string &index)
convert ChIndex name string to enum
PhiIndex
enum to classify the different phi layers in the muon spectrometer
constexpr int toInt(const EnumType enumVal)
ChIndex
enum to classify the different chamber layers in the muon spectrometer
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
Definition index.py:1
double deltaR2(double rapidity1, double phi1, double rapidity2, double phi2)
from bare rapidity,phi
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
Muon_v1 Muon
Reference the current persistent version:
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.