ATLAS Offline Software
Loading...
Searching...
No Matches
ActsMuonTrackingGeometryTest.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
9#include "GaudiKernel/EventContext.h"
10#include "GaudiKernel/ISvcLocator.h"
15
22
23
24
25//ACTS
26#include "Acts/Definitions/Units.hpp"
27#include "Acts/EventData/ParticleHypothesis.hpp"
28#include "Acts/Propagator/ActorList.hpp"
29#include "Acts/Propagator/Propagator.hpp"
30#include "Acts/Propagator/MaterialInteractor.hpp"
31#include "Acts/Propagator/EigenStepper.hpp"
32#include "Acts/Propagator/Navigator.hpp"
33#include "Acts/Propagator/detail/SteppingLogger.hpp"
34#include "Acts/Surfaces/StrawSurface.hpp"
35#include "Acts/Utilities/AngleHelpers.hpp"
36
37namespace {
38using SegLink_t = std::vector<ElementLink<xAOD::MuonSegmentContainer>>;
39static const SG::ConstAccessor<SegLink_t> segAcc{"truthSegmentLinks"};
40static const Amg::Vector3D dummyVec{100.*Gaudi::Units::m, 100.*Gaudi::Units::m, 100.*Gaudi::Units::m};
41
42struct PropagatorRecorder{
44 Amg::Vector3D actsPropPos{dummyVec};
46 Amg::Vector3D actsGlobalPos{dummyVec};
48 Amg::Vector3D actsPropDir{Amg::Vector3D::Zero()};
50 Identifier id{};
52 double actsPropabsMomentum{0.};
54 double actsStepSize{0.};
56 double actsHitWireDist{-1.};
57};
58
59}
60
61
62namespace ActsTrk {
63
64
66 const MuonGMR4::MuonReadoutElement* reElement = m_r4DetMgr->getReadoutElement(hitId);
67 const IdentifierHash trfHash = reElement->detectorType() == ActsTrk::DetectorType::Mdt ?
68 reElement->measurementHash(hitId) : reElement->layerHash(hitId);
69 return reElement->globalToLocalTransform(gctx, trfHash);
70 }
71
73 const MuonGMR4::MuonReadoutElement* reElement = m_r4DetMgr->getReadoutElement(hitId);
74 const IdentifierHash trfHash = reElement->detectorType() == ActsTrk::DetectorType::Mdt ?
75 reElement->measurementHash(hitId) : reElement->layerHash(hitId);
76 return reElement->localToGlobalTransform(gctx, trfHash);
77 }
78
79
81 return m_r4DetMgr->getReadoutElement(id)->layerHash(id);
82 }
83
84
86 ATH_CHECK(AthHistogramAlgorithm::initialize());
87 ATH_CHECK(m_idHelperSvc.retrieve());
89 ATH_CHECK(m_detMgrKey.initialize());
90 ATH_CHECK(m_rndmGenSvc.retrieve());
91 ATH_CHECK(m_geoCtxKey.initialize());
93 ATH_CHECK(m_truthParticleKey.initialize());
94 ATH_CHECK(m_truthSegLinkKey.initialize());
95 ATH_CHECK(detStore()->retrieve(m_r4DetMgr));
96
97 ATH_CHECK(m_tree.init(this));
98 ATH_MSG_INFO("ActsMuonTrackingGeometryTest successfully initialized");
99 return StatusCode::SUCCESS;
100 }
101
103 ATH_CHECK(m_tree.write());
104 return StatusCode::SUCCESS;
105 }
106
108
109 const EventContext& ctx = Gaudi::Hive::currentContext();
110
111 const ActsTrk::GeometryContext* gctx{nullptr};
112 const AtlasFieldCacheCondObj* fieldCondObj{nullptr};
113 const MuonGM::MuonDetectorManager* detMgr{nullptr};
114 const xAOD::TruthParticleContainer* truthParticles{nullptr};
115
116 ATH_CHECK(SG::get(fieldCondObj, m_fieldCacheCondObjInputKey, ctx));
117 ATH_CHECK(SG::get(detMgr, m_detMgrKey, ctx));
118 ATH_CHECK(SG::get(truthParticles, m_truthParticleKey, ctx));
119 ATH_CHECK(SG::get(gctx, m_geoCtxKey, ctx));
120
121
122 const Acts::MagneticFieldContext mfContext = Acts::MagneticFieldContext(fieldCondObj);
123
124 auto anygctx = gctx->context();
125
126 //Get the tracking geometry
127 auto trackingGeometry = m_trackingGeometryTool->trackingGeometry();
128
129 if (!trackingGeometry) {
130 ATH_MSG_ERROR("Failed to retrieve the tracking geometry");
131 return StatusCode::FAILURE;
132 }
133
134 //Configure the ACTS propagator with the navigator and the stepper
135 using Stepper = Acts::EigenStepper<>;
136 using Navigator = Acts::Navigator;
137 using Propagator = Acts::Propagator<Stepper,Navigator>;
138 using ActorList = Acts::ActorList<Acts::detail::SteppingLogger, Acts::MaterialInteractor, Acts::EndOfWorldReached>;
139 using PropagatorOptions = Propagator::Options<ActorList>;
140
141 Navigator::Config navCfg;
142 navCfg.trackingGeometry = trackingGeometry;
143 navCfg.resolveSensitive = true;
144 navCfg.resolveMaterial= false;
145 navCfg.resolvePassive = true;
146
147 Navigator navigator(navCfg, Acts::getDefaultLogger("Navigator", ActsTrk::actsLevelVector(msgLevel())));
148
149 auto bField = std::make_shared<ATLASMagneticFieldWrapper>();
150 auto bfieldCache = bField->makeCache(mfContext);
151
152 auto stepper = Stepper(bField);
153
154 PropagatorOptions options(anygctx, mfContext);
155 options.pathLimit = m_pathLimit;
156 options.stepping.maxStepSize = m_maxStepSize;
157 options.stepping.stepTolerance = m_stepTolerance;
158 options.maxSteps = m_maxSteps;
159 options.maxTargetSkipping = m_maxTargetSkipping;
160
161
162 //switch off material interactions
163 auto& materialInteractor = options.actorList.get<Acts::MaterialInteractor>();
164 materialInteractor.energyLoss = false;
165 materialInteractor.multipleScattering = false;
166 materialInteractor.recordInteractions = false;
167
168 Propagator propagator(std::move(stepper), std::move(navigator),
169 Acts::getDefaultLogger("Propagator", ActsTrk::actsLevelVector(msgLevel())));
170
171 for(const xAOD::TruthParticle* truthParticle : *truthParticles) {
172
173 ATH_MSG_DEBUG("Processing truth particle with PDG ID: " << truthParticle->pdgId() << " ,pT: "
174 << truthParticle->pt() << " , p: " << truthParticle->p4().P() << ", eta: " << truthParticle->eta() << " , phi: " << truthParticle->phi());
175
176 //select only muons or geantinos particles
177 if(std::abs(truthParticle->pdgId()) != 13 && std::abs(truthParticle->pdgId()) != 998) {
178 ATH_MSG_VERBOSE("Skipping truth particle with PDG ID: " << truthParticle->pdgId()<<" only muons or charged geantinos are being processed");
179 continue;
180 }
181
182 Acts::ParticleHypothesis actsParticleHypothesis = truthParticle->pdgId() == 998 ?
183 Acts::ParticleHypothesis::chargedGeantino() : Acts::ParticleHypothesis::muon();
184
185
186 //take the truth segments witht he sim hits
187 std::vector<std::pair<const xAOD::MuonSegment*, std::vector<const xAOD::MuonSimHit*>>> muonSegmentWithSimHits;
188 const SegLink_t& segLink = segAcc(*truthParticle);
189 if (segLink.empty()) {
190 ATH_MSG_WARNING("No segment link found for truth particle with PDG ID: " << truthParticle->pdgId());
191 continue;
192 }
193 for(const auto& truthSegLink: segLink) {
194 const xAOD::MuonSegment* seg{*truthSegLink};
195 if(!seg){
196 continue;
197 }
198 auto unordedHits = MuonR4::getMatchingSimHits(*seg);
199 std::vector<const xAOD::MuonSimHit*> muonSimHits{unordedHits.begin(), unordedHits.end()};
200 ATH_MSG_VERBOSE("Segment at : "<<Amg::toString(seg->position())<<" with Sim Hits: "<<unordedHits.size());
201 muonSegmentWithSimHits.emplace_back(seg,muonSimHits);
202
203 }
204
205 const auto& particle = truthParticle->p4();
206 m_truthPt = particle.Pt();
207 m_truthP = particle.P();
208
209 const xAOD::TruthVertex* prodVertex = truthParticle->prodVtx();
210 Amg::Vector3D prodPos = prodVertex ? Amg::Vector3D(prodVertex->x(), prodVertex->y(), prodVertex->z())
211 : Amg::Vector3D::Zero();
212 Amg::Vector3D prodDir{Amg::dirFromAngles(particle.Phi(), particle.Theta())};
213 ATH_MSG_VERBOSE("Truth particle produced at " << Amg::toString(prodPos) << " with direction " << Amg::toString(prodDir));
214
215 //position, direction and momentum at the start of the propagation
216 Amg::Vector3D startPropPos{prodPos};
217 Amg::Vector3D startPropDir{prodDir};
218 double startPropP{energyToActs(particle.P())};
219
221
222 if(muonSegmentWithSimHits.empty()) {
223 ATH_MSG_DEBUG("No segments found for truth particle with PDG ID: " << truthParticle->pdgId());
224 continue;
225 }
226 // sort the segments by their radial distance
227 std::ranges::sort(muonSegmentWithSimHits,
228 [](const auto& a, const auto& b) {
229 return a.first->position().perp() < b.first->position().perp();
230 });
231
232
233 //sort the sim hits of the first segment to get the starting propagation position as the first one
234 std::vector<const xAOD::MuonSimHit*>& muonSimHits = muonSegmentWithSimHits.front().second;
235 std::ranges::sort(muonSimHits,
236 [this, &gctx](const xAOD::MuonSimHit* hit1, const xAOD::MuonSimHit* hit2){
237 const Amg::Vector3D globalPos1 = toGlobalTrf(*gctx,hit1->identify())*xAOD::toEigen(hit1->localPosition());
238 const Amg::Vector3D globalPos2 = toGlobalTrf(*gctx,hit2->identify())*xAOD::toEigen(hit2->localPosition());
239 return globalPos1.norm()<globalPos2.norm();
240 });
241
242 ATH_MSG_VERBOSE("After sorting..");
243 for(const auto& simHit: muonSimHits){
244 ATH_MSG_VERBOSE("Sim Hit of first segment at position: "<<Amg::toString(toGlobalTrf(*gctx, simHit->identify())*xAOD::toEigen(simHit->localPosition())));
245
246 }
247
248
249 startPropPos = toGlobalTrf(*gctx, muonSimHits.front()->identify())*xAOD::toEigen(muonSimHits.front()->localPosition());
250 startPropDir = toGlobalTrf(*gctx, muonSimHits.front()->identify()).linear()*xAOD::toEigen(muonSimHits.front()->localDirection());
251 startPropP = energyToActs(muonSimHits.front()->kineticEnergy());
252
253
254 ATH_MSG_VERBOSE("Kinetic Energy from the simHit: "<<startPropP / Gaudi::Units::GeV<<" and mass: "<<muonSimHits.front()->mass()<<" and energy deposit: "<<muonSimHits.front()->energyDeposit() / Gaudi::Units::eV);
255
256
257 }
258
259 ATH_MSG_VERBOSE("Starting propagation from"<<Amg::toString(startPropPos)<<" with direction"
260 <<Amg::toString(startPropDir)<<" and momentum "<<startPropP);
261
262
263
264 Acts::BoundTrackParameters start = Acts::BoundTrackParameters::createCurvilinear(
265 Acts::VectorHelpers::makeVector4(startPropPos, 0.), startPropDir,
266 truthParticle->charge() / startPropP,
267 std::nullopt,
268 actsParticleHypothesis);
269
270 ATH_MSG_DEBUG("start propagating here");
271 auto result = propagator.propagate(start, options);
272 const Acts::detail::SteppingLogger::result_type state = result.value().get<Acts::detail::SteppingLogger::result_type>();
273 const Acts::MaterialInteractor::result_type material = result.value().get<Acts::MaterialInteractor::result_type>();
274
275 m_propSteps = state.steps.size();
276 ATH_MSG_DEBUG("Number of propagated steps : " << m_propSteps);
277 m_propLength = result.value().pathLength;
278 std::vector<PropagatorRecorder> propagatedHits;
279
280 for(const auto& step : state.steps) {
281 if(!step.surface) {
282 continue;
283 }
284
285 const SurfaceCache* sCache = dynamic_cast<const SurfaceCache *>(step.surface->associatedDetectorElement());
286 if(!sCache) {
287 ATH_MSG_VERBOSE("Surface found but it's a portal, continuing..");
288 continue;
289 }
290 const Identifier ID = sCache->identify();
291 const Amg::Transform3D toGap{toLocalTrf(*gctx, ID)};
292 ATH_MSG_VERBOSE("Identify propagated hit " << m_idHelperSvc->toString(ID) << " with hit at local position " << Amg::toString(toGap*step.position)<<" and global direction "<<Amg::toString(step.momentum.unit()));
293
294
295 PropagatorRecorder newRecord{};
296 newRecord.id = ID;
297 newRecord.actsPropPos = toGap*step.position;
298 newRecord.actsGlobalPos = step.position;
299 newRecord.actsPropDir = toGap.linear()*step.momentum.unit();
300 newRecord.actsPropabsMomentum = step.momentum.norm();
301 newRecord.actsStepSize = step.stepSize.value();
302 //calculate the distance of the propagated hit to the wire in case of MDTs
303
304 if(m_idHelperSvc->isMdt(ID)){
305 //get the surface
306 const auto* tubeSurf = dynamic_cast<const Acts::StrawSurface*>(&sCache->surface());
307 if(tubeSurf){
308 Amg::Vector3D wirePos = toGap*tubeSurf->center(anygctx);
309 Amg::Vector3D wireDir = toGap.linear()*tubeSurf->lineDirection(anygctx);
310 Amg::Vector3D localPropHit = newRecord.actsPropPos;
311 Amg::Vector3D dirPropHit = newRecord.actsPropDir;
312 double distToWire = Amg::lineDistance(wirePos, wireDir, localPropHit, dirPropHit);
313 newRecord.actsHitWireDist = distToWire;
314 }
315
316 }
317
318 propagatedHits.emplace_back(newRecord);
319
320 }
321
322
323 for (const auto&[segment, simHits] : muonSegmentWithSimHits) {
324 for(const xAOD::MuonSimHit* simHit : simHits){
325
326 const Identifier ID = simHit->identify();
327
328 const Amg::Vector3D localPos = xAOD::toEigen(simHit->localPosition());
329
330 const Amg::Vector3D globalPos = toGlobalTrf(*gctx, simHit->identify())*xAOD::toEigen(simHit->localPosition());
331 const Amg::Vector3D localDir = xAOD::toEigen(simHit->localDirection());
332
333 if(m_r4DetMgr->getReadoutElement(ID)->detectorType() == ActsTrk::DetectorType::sTgc && localPos.z() != 0.0){
334 continue;
335
336 }
337
338 m_detId.push_back(ID);
339 m_techIdx.push_back(toInt(m_idHelperSvc->technologyIndex(ID)));
340 m_gasGapId.push_back(layerHash(ID));
341
342 m_truthLoc.push_back(localPos);
343 m_truthDir.push_back(localDir);
344 m_truthGlob.push_back(globalPos);
345 m_startGlob.push_back(startPropPos);
346 ATH_MSG_DEBUG("Truth hit with ID: " << m_idHelperSvc->toString(ID)
347 << " at local position: " << Amg::toString(localPos)
348 << " and global position: " << Amg::toString(globalPos)
349 << " and direction: " << Amg::toString(localDir));
350
351 //get the matching propagated hits
352 auto it_begin = std::ranges::find_if(propagatedHits,
353 [this, ID](const auto& propagatedHit) {
354 return m_idHelperSvc->detElId(ID) == m_idHelperSvc->detElId(propagatedHit.id) &&
355 layerHash(ID) == layerHash(propagatedHit.id);
356 });
357 m_isPropagated.push_back(it_begin != propagatedHits.end());
358
359 if(it_begin == propagatedHits.end()){
360 m_actsPropLoc.push_back(dummyVec);
361 Amg::Vector3D zero = Amg::Vector3D::Zero();
362 m_actsPropDir.push_back(zero);
363 m_actsPropGlob.push_back(dummyVec);
364 m_actsPropabsMomentum.push_back(0.);
365 m_actsStepSize.push_back(0.);
366 m_actsHitWireDist.push_back(-1.);
367 continue;
368 }
369
370 //find - if any - propagated hits on the same layer
371 auto it_end = std::find_if(it_begin, propagatedHits.end(),
372 [this, ID](const auto& propagatedHit) {
373 return m_idHelperSvc->detElId(ID) != m_idHelperSvc->detElId(propagatedHit.id) ||
374 layerHash(ID) != layerHash(propagatedHit.id);
375 });
376
378 auto it = std::min_element(it_begin, it_end,
379 [&localPos](const PropagatorRecorder& a,
380 const PropagatorRecorder& b){
381 return (localPos - a.actsPropPos).mag() < (localPos - b.actsPropPos).mag();
382 });
383
384 ATH_MSG_DEBUG("Found propagated hit with ID: " << m_idHelperSvc->toString(it->id)
385 << " at local position: " << Amg::toString(it->actsPropPos)
386 << " and global position: " << Amg::toString(it->actsGlobalPos)
387 << " and direction: " << Amg::toString(it->actsPropDir));
388
389 m_actsPropLoc.push_back(it->actsPropPos);
390 m_actsPropDir.push_back(it->actsPropDir);
391 m_actsPropGlob.push_back(it->actsGlobalPos);
392 m_actsPropabsMomentum.push_back(it->actsPropabsMomentum);
393 m_actsStepSize.push_back(it->actsStepSize);
394 m_actsHitWireDist.push_back(it->actsHitWireDist);
395
396
397 }
398 }
399
400 m_event = ctx.eventID().event_number();
401 m_tree.fill(ctx);
402
403
404 }
405
406
407 return StatusCode::SUCCESS;
408 }
409
410}
411
Scalar mag() const
mag method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::vector< Identifier > ID
ATLAS-specific HepMC functions.
static Double_t a
MuonVal::ScalarBranch< unsigned int > & m_propSteps
MuonVal::VectorBranch< unsigned short > & m_isPropagated
SG::ReadHandleKey< ActsTrk::GeometryContext > m_geoCtxKey
MuonVal::VectorBranch< unsigned short > & m_gasGapId
MuonVal::VectorBranch< float > & m_actsHitWireDist
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthParticleKey
const MuonGMR4::MuonDetectorManager * m_r4DetMgr
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
Amg::Transform3D toGlobalTrf(const ActsTrk::GeometryContext &gctx, const Identifier &hitId) const
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_truthSegLinkKey
MuonVal::VectorBranch< float > & m_actsPropabsMomentum
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
IdentifierHash layerHash(const Identifier &id) const
MuonVal::ScalarBranch< unsigned int > & m_event
Amg::Transform3D toLocalTrf(const ActsTrk::GeometryContext &gctx, const Identifier &hitId) const
MuonVal::VectorBranch< unsigned short > & m_techIdx
SG::ReadCondHandleKey< MuonGM::MuonDetectorManager > m_detMgrKey
MuonVal::VectorBranch< float > & m_actsStepSize
PublicToolHandle< ActsTrk::ITrackingGeometryTool > m_trackingGeometryTool
Acts::GeometryContext context() const
virtual DetectorType detectorType() const =0
Returns the detector element type.
: Helper class to connect the aligned transformations of each active sensor(layer) with the Acts::Sur...
Identifier identify() const override final
Returns the identifier of the Surface.
const Acts::Surface & surface() const override final
Returns the dereferenced pointer cache.
const ServiceHandle< StoreGateSvc > & detStore() const
This is a "hash" representation of an Identifier.
MuonReadoutElement is an abstract class representing the geometry of a muon detector.
const Amg::Transform3D & localToGlobalTransform(const ActsTrk::GeometryContext &ctx) const
Returns the transformation from the local coordinate system of the readout element into the global AT...
Amg::Transform3D globalToLocalTransform(const ActsTrk::GeometryContext &ctx) const
Returns the transformation from the global ATLAS coordinate system into the local coordinate system o...
virtual IdentifierHash measurementHash(const Identifier &measId) const =0
The measurement hash is a continous numbering schema of all readout channels described by the specifi...
virtual IdentifierHash layerHash(const Identifier &measId) const =0
The layer hash removes the bits from the IdentifierHash corresponding to the measurement's channel nu...
The MuonDetectorManager stores the transient representation of the Muon Spectrometer geometry and pro...
Helper class to provide constant type-safe access to aux data.
Amg::Vector3D position() const
Returns the position as Amg::Vector.
Identifier identify() const
Returns the global ATLAS identifier of the SimHit.
ConstVectorMap< 3 > localPosition() const
Returns the local postion of the traversing particle.
float z() const
Vertex longitudinal distance along the beam line form the origin.
float y() const
Vertex y displacement.
float x() const
Vertex x displacement.
void zero(TH2 *h)
zero the contents of a 2d histogram
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
constexpr double energyToActs(const double athenaE)
Converts an energy scalar from Athena to Acts units.
Acts::Logging::Level actsLevelVector(MSG::Level lvl)
@ sTgc
Micromegas (NSW)
@ Mdt
MuonSpectrometer.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
double lineDistance(const AmgVector(N)&posA, const AmgVector(N)&dirA, const AmgVector(N)&posB, const AmgVector(N)&dirB)
: Calculates the shortest distance between two lines
Amg::Vector3D dirFromAngles(const double phi, const double theta)
Constructs a direction vector from the azimuthal & polar angles.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
std::unordered_set< const xAOD::MuonSimHit * > getMatchingSimHits(const xAOD::MuonSegment &segment)
: Returns all sim hits matched to a xAOD::MuonSegment
ElementLink< MuonR4::SegmentContainer > SegLink_t
Abrivation of the link to the reco segment container.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
MuonSimHit_v1 MuonSimHit
Defined the version of the MuonSimHit.
Definition MuonSimHit.h:12
TruthParticle_v1 TruthParticle
Typedef to implementation.
MuonSegment_v1 MuonSegment
Reference the current persistent version:
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.