ATLAS Offline Software
Loading...
Searching...
No Matches
MuonSensitiveDetector.cxx
Go to the documentation of this file.
1
2/*
3 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
4*/
5
7
10#include <GeoModelKernel/throwExcept.h>
12
13#include <G4Geantino.hh>
14#include <G4ChargedGeantino.hh>
15
16#include <MCTruth/TrackHelper.h>
19
20using namespace ActsTrk;
21
22namespace {
23 static const SG::Decorator<int> dec_G4TrkId{"MuonSim_G4TrkId"};
24}
25
26namespace MuonG4R4 {
28 const std::string& output_key,
29 const std::string& trfStore_key,
30 const MuonGMR4::MuonDetectorManager* detMgr):
31 G4VSensitiveDetector{name},
32 AthMessaging{name},
33 m_writeKey{output_key},
34 m_trfCacheKey{trfStore_key},
35 m_detMgr{detMgr} {
36 m_trfCacheKey.initialize().ignore();
37 }
38 void MuonSensitiveDetector::Initialize(G4HCofThisEvent*) {
41 auto* hitVec = m_g4UserEventInfo->GetHitCollectionMap()->Find<MuonSimHitsVec>(m_writeKey);
42 if (!hitVec) {
43 THROW_EXCEPTION("The event does not contain a MuonSimHit container called '"<<m_writeKey<<"'");
44 }
45 m_outContainer = hitVec->container.get();
46 ATH_MSG_DEBUG(__func__<<"() "<<__LINE__<<" - Retrieved '"<<m_writeKey<<"'.");
47 } else {
48 m_outContainer = nullptr;
49 }
50 }
54 const EventContext& MuonSensitiveDetector::eventContext() const {
55 AtlasG4EventUserInfo* g4EventInfo = eventInfo();
56 if (!g4EventInfo) {
57 THROW_EXCEPTION("No AtlasG4EventUserInfo available");
58 }
59 return g4EventInfo->GetEventContext();
60 }
63 const ActsTrk::DetectorAlignStore* alignment{nullptr};
64 if (!SG::get(alignment, m_trfCacheKey, eventContext()).isSuccess()) {
65 THROW_EXCEPTION("Failed to retrieve "<<m_trfCacheKey.fullKey()<<".");
66 }
67 gctx.setStore(std::make_unique<DetectorAlignStore>(*alignment));
68 return gctx;
69 }
70 bool MuonSensitiveDetector::processStep(const G4Step* aStep) const {
71 const G4Track* currentTrack = aStep->GetTrack();
72 ATH_MSG_VERBOSE("Check whether step pdgId: "<<(*currentTrack)<<" will be processed.");
74 constexpr double velCutOff = 10.*Gaudi::Units::micrometer / Gaudi::Units::second;
75 if (aStep->GetStepLength() < std::numeric_limits<float>::epsilon() || currentTrack->GetVelocity() < velCutOff) {
76 ATH_MSG_VERBOSE("Step length is too short ");
77 return false;
78 }
80 if (currentTrack->GetDefinition()->GetPDGCharge() == 0.0) {
81 ATH_MSG_VERBOSE("Particle is neutral");
82 return currentTrack->GetDefinition() == G4Geantino::GeantinoDefinition() ||
83 currentTrack->GetDefinition() == G4ChargedGeantino::ChargedGeantinoDefinition();
84 }
85 return true;
86 }
88 TrackHelper trHelper{track};
89 AtlasG4EventUserInfo* g4EventInfo = eventInfo();
90 return trHelper.GenerateParticleLink(g4EventInfo ? g4EventInfo->GetEventStore() : nullptr);
91 }
93 const Amg::Transform3D& toGasGap,
94 const G4Step* aStep) {
95
96 const G4Track* currentTrack = aStep->GetTrack();
98 const Amg::Vector3D locPostStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPostStepPoint()->GetPosition())};
99 Amg::Vector3D locPreStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPreStepPoint()->GetPosition())};
100
102 Amg::Vector3D locHitDir = toGasGap.linear() * Amg::Hep3VectorToEigen(currentTrack->GetMomentumDirection());
103
106 xAOD::MuonSimHit* prevHit = lastSnapShot(hitID, aStep);
107 if (std::abs(currentTrack->GetParticleDefinition()->GetPDGEncoding()) == 11 && prevHit) {
108 locPreStep = xAOD::toEigen(prevHit->localPosition()) - 0.5* prevHit->stepLength() *
109 xAOD::toEigen(prevHit->localDirection());
110
111 locHitDir = (locPostStep - locPreStep).unit();
112 }
113 const Amg::Vector3D locHitPos = 0.5* (locPreStep + locPostStep);
114 ATH_MSG_VERBOSE( m_detMgr->idHelperSvc()->toStringGasGap(hitID)<<" - track: "<<(*currentTrack)
115 <<", deposit: "<<aStep->GetTotalEnergyDeposit()<<", -- local coords: "
116 <<"prestep: "<<Amg::toString(locPreStep)<<", post step: "<<Amg::toString(locPostStep)
117 <<" mid point: "<< Amg::toString(locHitPos)<<", direction: "<<Amg::toString(locHitDir)
118 <<", deposit: "<<aStep->GetTotalEnergyDeposit());
119
120 const double globalTime = currentTrack->GetGlobalTime() + locHitDir.dot(locPostStep - locHitPos) / currentTrack->GetVelocity();
121
122 xAOD::MuonSimHit* newHit = saveHit(hitID, locHitPos, locHitDir, globalTime, aStep);
123 if (prevHit) {
124 newHit->setStepLength((locPostStep - locPreStep).mag());
125 }
126 return newHit;
127 }
129 const Amg::Vector3D& hitPos,
130 const Amg::Vector3D& hitDir,
131 const double globTime,
132 const G4Step* aStep) {
133 const G4Track* currentTrack = aStep->GetTrack();
134 const HepMcParticleLink particleLink{genParticleLink(currentTrack)};
135 // If Geant4 propagates the same track through the same volume just update the last hit and don't write a new one
136 xAOD::MuonSimHit* hit = lastSnapShot(hitId, aStep, particleLink);
137 bool newHit{false};
138 if (!hit) {
139 if (!m_outContainer) {
140 THROW_EXCEPTION("Cannot find container '"<<m_writeKey<<"'.");
141 }
142 hit = m_outContainer->push_back(std::make_unique<xAOD::MuonSimHit>());
143 newHit = true;
144 }
145 dec_G4TrkId(*hit) = currentTrack->GetTrackID();
146 hit->setIdentifier(hitId);
147 hit->setLocalPosition(xAOD::toStorage(hitPos));
149 hit->setMass(currentTrack->GetDefinition()->GetPDGMass());
150 hit->setGlobalTime(globTime);
151 hit->setPdgId(currentTrack->GetDefinition()->GetPDGEncoding());
152 hit->setEnergyDeposit(aStep->GetTotalEnergyDeposit() + (newHit ? 0. : hit->energyDeposit()));
153 hit->setKineticEnergy(currentTrack->GetKineticEnergy());
154 hit->setGenParticleLink(particleLink);
155 hit->setStepLength(aStep->GetStepLength());
156
157 ATH_MSG_VERBOSE("Save new hit "<<m_detMgr->idHelperSvc()->toString(hitId)
158 <<", pdgId: "<<hit->pdgId()
159 <<", "<<particleLink
160 <<", trackId: "<<currentTrack->GetTrackID()<<", "
161 <<", "<<particleLink.cptr()<<std::endl
162 <<"pos: "<<Amg::toString(hitPos)<<", dir: "<<Amg::toString(hitDir)<<", time: "<<globTime
163 <<", energy: "<<hit->kineticEnergy()<<", stepLength: "<<hit->stepLength()<<", "
164 <<", deposit energy: "<<hit->energyDeposit());
165 return hit;
166 }
168 const G4Step* hitStep) {
169 const HepMcParticleLink particleLink{genParticleLink(hitStep->GetTrack())};
170 return lastSnapShot(hitId, hitStep, particleLink);
171 }
173 const G4Step* hitStep,
174 const HepMcParticleLink& particleLink) {
177 if (m_outContainer->empty() ||
178 m_outContainer->back()->identify() != hitId ||
179 particleLink != m_outContainer->back()->genParticleLink() ||
180 dec_G4TrkId(*m_outContainer->back()) != hitStep->GetTrack()->GetTrackID()) {
181 return nullptr;
182 }
183 return m_outContainer->back();
184 }
185}
Scalar mag() const
mag method
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Utility helpers used by Run-4 muon sensitive detector implementations.
void setStore(AlignmentStorePtr store)
Adds the store to the Geometry context.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
This class is attached to G4Event objects as UserInformation.
static AtlasG4EventUserInfo * GetEventUserInfo()
const EventContext & GetEventContext() const
const T * get(size_type n) const
Access an element, as an rvalue.
const MuonGMR4::MuonDetectorManager * m_detMgr
Pointer to the underlying detector manager.
AtlasG4EventUserInfo * eventInfo() const
Returns the event user info for the current G4 event, if available.
const EventContext & eventContext() const
Returns the current event context recorded in the G4 event info.
virtual void Initialize(G4HCofThisEvent *HCE) override final
Create the output container at the beginning of the event.
MuonSensitiveDetector(const std::string &name, const std::string &output_key, const std::string &trf_storeKey, const MuonGMR4::MuonDetectorManager *detMgr)
Constructor.
SG::ReadHandleKey< ActsTrk::DetectorAlignStore > m_trfCacheKey
ReadHandleKey to the DetectorAlignmentStore caching the relevant transformations needed in this event...
HepMcParticleLink genParticleLink(const G4Track *track) const
Generates a HepMcParticleLink for the Geant4 track using the current event store.
xAOD::MuonSimHit * propagateAndSaveStrip(const Identifier &hitId, const Amg::Transform3D &toGasGap, const G4Step *hitStep)
Records the G4Step in the sim hit.
bool processStep(const G4Step *step) const
Checks whether the current step shall be processed at all.
AtlasG4EventUserInfo * m_g4UserEventInfo
G4 event user information holding the current event store/context.
std::string m_writeKey
Key under which the output container is stored in the G4 event.
xAOD::MuonSimHit * saveHit(const Identifier &hitId, const Amg::Vector3D &hitPos, const Amg::Vector3D &hitDir, const double globTime, const G4Step *hitStep)
Saves the current Step as a xAOD::MuonSimHit snapshot.
ActsTrk::GeometryContext getGeoContext() const
Returns the current geometry context in the event.
xAOD::MuonSimHitContainer * m_outContainer
Pointer to the MuonSimHit output container.
xAOD::MuonSimHit * lastSnapShot(const Identifier &gasGapId, const G4Step *hitStep, const HepMcParticleLink &particleLink)
Returns the last snap shot matching an already generated HepMcParticleLink.
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
void setGenParticleLink(const HepMcParticleLink &link)
Sets the link to the HepMC particle producing this hit.
void setGlobalTime(const float time)
Sets the time of the traversing particle.
void setLocalDirection(MeasVector< 3 > vec)
Sets the local direction of the traversing particle.
void setEnergyDeposit(const float deposit)
Sets the energy deposited by the traversing particle inside the gas volume.
void setMass(const float m)
set the rest-mass of the traversing particle
ConstVectorMap< 3 > localDirection() const
Returns the local direction of the traversing particle.
int pdgId() const
Returns the pdgID of the traversing particle.
void setIdentifier(const Identifier &id)
Sets the global ATLAS identifier.
float stepLength() const
Returns the path length of the G4 step.
void setStepLength(const float length)
Set the path length of the G4 step.
float energyDeposit() const
Returns the energy deposited by the traversing particle inside the gas volume.
ConstVectorMap< 3 > localPosition() const
Returns the local postion of the traversing particle.
void setKineticEnergy(const float energy)
Sets the kinetic energy of the traversing particle.
void setLocalPosition(MeasVector< 3 > vec)
Sets the local position of the traversing particle.
void setPdgId(int id)
Sets the pdgID of the traversing particle.
float kineticEnergy() const
Returns the kinetic energy of the traversing particle.
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Amg::Vector3D Hep3VectorToEigen(const CLHEP::Hep3Vector &CLHEPvector)
Converts a CLHEP-based CLHEP::Hep3Vector into an Eigen-based Amg::Vector3D.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Include the common definitions from the MuonReadoutGeometry.
AthHitVec::AuxStoreHitCollection< xAOD::MuonSimHitContainer, xAOD::MuonSimHitAuxContainer > MuonSimHitsVec
Helper type which can be filled into the HitCollectionMap storing the event content in the event.
SG::Decorator< T, ALLOC > Decorator
Helper class to provide type-safe access to aux data, specialized for JaggedVecElt.
Definition AuxElement.h:576
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
MuonSimHit_v1 MuonSimHit
Defined the version of the MuonSimHit.
Definition MuonSimHit.h:12
MeasVector< N > toStorage(const AmgVector(N)&amgVec)
Converts the double precision of the AmgVector into the floating point storage precision of the MeasV...
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10