ATLAS Offline Software
Loading...
Searching...
No Matches
MuonSensitiveDetector.cxx
Go to the documentation of this file.
1
2/*
3 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
4*/
5
7
10#include <GeoModelKernel/throwExcept.h>
12
13#include <MCTruth/TrackHelper.h>
14#include <G4Geantino.hh>
15#include <G4ChargedGeantino.hh>
16
18
19using namespace ActsTrk;
20
21namespace {
22 static const SG::Decorator<int> dec_G4TrkId{"MuonSim_G4TrkId"};
23}
24
25namespace MuonG4R4 {
27 const std::string& output_key,
28 const std::string& trfStore_key,
30 G4VSensitiveDetector{name},
31 AthMessaging{name},
32 m_writeHandle{output_key},
33 m_trfCacheKey{trfStore_key},
34 m_detMgr{detMgr} {
35 m_trfCacheKey.initialize().ignore();
36 }
37 void MuonSensitiveDetector::Initialize(G4HCofThisEvent*) {
38 if (m_writeHandle.isValid()) {
39 ATH_MSG_VERBOSE("Simulation hit container "<<m_writeHandle.fullKey()<<" is already written");
40 return;
41 }
42 if (!m_writeHandle.recordNonConst(std::make_unique<xAOD::MuonSimHitContainer>(),
43 std::make_unique<xAOD::MuonSimHitAuxContainer>()).isSuccess()) {
44 THROW_EXCEPTION(" Failed to record "<<m_writeHandle.fullKey());
45 }
46 ATH_MSG_DEBUG("Output container "<<m_writeHandle.fullKey()<<" has been successfully created");
47 }
50 SG::ReadHandle trfStoreHandle{m_trfCacheKey};
51 if (!trfStoreHandle.isValid()) {
52 THROW_EXCEPTION("Failed to retrieve "<<m_trfCacheKey.fullKey()<<".");
53 }
54 gctx.setStore(std::make_unique<DetectorAlignStore>(*trfStoreHandle));
55 return gctx;
56 }
57 bool MuonSensitiveDetector::processStep(const G4Step* aStep) const {
58 const G4Track* currentTrack = aStep->GetTrack();
59 ATH_MSG_VERBOSE("Check whether step pdgId: "<<(*currentTrack)<<" will be processed.");
61 constexpr double velCutOff = 10.*Gaudi::Units::micrometer / Gaudi::Units::second;
62 if (aStep->GetStepLength() < std::numeric_limits<float>::epsilon() || currentTrack->GetVelocity() < velCutOff) {
63 ATH_MSG_VERBOSE("Step length is too short ");
64 return false;
65 }
67 if (currentTrack->GetDefinition()->GetPDGCharge() == 0.0) {
68 ATH_MSG_VERBOSE("Particle is neutral");
69 return currentTrack->GetDefinition() == G4Geantino::GeantinoDefinition() ||
70 currentTrack->GetDefinition() == G4ChargedGeantino::ChargedGeantinoDefinition();
71 }
72 return true;
73 }
75 const Amg::Transform3D& toGasGap,
76 const G4Step* aStep) {
77
78 const G4Track* currentTrack = aStep->GetTrack();
80 const Amg::Vector3D locPostStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPostStepPoint()->GetPosition())};
81 Amg::Vector3D locPreStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPreStepPoint()->GetPosition())};
82
84 Amg::Vector3D locHitDir = toGasGap.linear() * Amg::Hep3VectorToEigen(currentTrack->GetMomentumDirection());
85
88 xAOD::MuonSimHit* prevHit = lastSnapShot(hitID, aStep);
89 if (std::abs(currentTrack->GetParticleDefinition()->GetPDGEncoding()) == 11 && prevHit) {
90 locPreStep = xAOD::toEigen(prevHit->localPosition()) - 0.5* prevHit->stepLength() *
91 xAOD::toEigen(prevHit->localDirection());
92
93 locHitDir = (locPostStep - locPreStep).unit();
94 }
95 const Amg::Vector3D locHitPos = 0.5* (locPreStep + locPostStep);
96 ATH_MSG_VERBOSE( m_detMgr->idHelperSvc()->toStringGasGap(hitID)<<" - track: "<<(*currentTrack)
97 <<", deposit: "<<aStep->GetTotalEnergyDeposit()<<", -- local coords: "
98 <<"prestep: "<<Amg::toString(locPreStep)<<", post step: "<<Amg::toString(locPostStep)
99 <<" mid point: "<< Amg::toString(locHitPos)<<", direction: "<<Amg::toString(locHitDir)
100 <<", deposit: "<<aStep->GetTotalEnergyDeposit());
101
102 const double globalTime = currentTrack->GetGlobalTime() + locHitDir.dot(locPostStep - locHitPos) / currentTrack->GetVelocity();
103
104 xAOD::MuonSimHit* newHit = saveHit(hitID, locHitPos, locHitDir, globalTime, aStep);
105 if (prevHit) {
106 newHit->setStepLength((locPostStep - locPreStep).mag());
107 }
108 return newHit;
109 }
111 const Amg::Vector3D& hitPos,
112 const Amg::Vector3D& hitDir,
113 const double globTime,
114 const G4Step* aStep) {
115 const G4Track* currentTrack = aStep->GetTrack();
116 TrackHelper trHelper{currentTrack};
117 // If Geant4 propagates the same track through the same volume just update the last hit and don't write a new one
118 xAOD::MuonSimHit* hit = lastSnapShot(hitId, aStep);
119 bool newHit{false};
120 if (!hit) {
121 hit = m_writeHandle->push_back(std::make_unique<xAOD::MuonSimHit>());
122 newHit = true;
123 }
124 dec_G4TrkId(*hit) = currentTrack->GetTrackID();
125 hit->setIdentifier(hitId);
126 hit->setLocalPosition(xAOD::toStorage(hitPos));
128 hit->setMass(currentTrack->GetDefinition()->GetPDGMass());
129 hit->setGlobalTime(globTime);
130 hit->setPdgId(currentTrack->GetDefinition()->GetPDGEncoding());
131 hit->setEnergyDeposit(aStep->GetTotalEnergyDeposit() + (newHit ? 0. : hit->energyDeposit()));
132 hit->setKineticEnergy(currentTrack->GetKineticEnergy());
134 hit->setStepLength(aStep->GetStepLength());
135
136 ATH_MSG_VERBOSE("Save new hit "<<m_detMgr->idHelperSvc()->toString(hitId)
137 <<", pdgId: "<<hit->pdgId()
138 <<", "<<hit->genParticleLink()
139 <<", trackId: "<<currentTrack->GetTrackID()<<", "
140 <<", "<<hit->genParticleLink().cptr()<<std::endl
141 <<"pos: "<<Amg::toString(hitPos)<<", dir: "<<Amg::toString(hitDir)<<", time: "<<globTime
142 <<", energy: "<<hit->kineticEnergy()<<", stepLength: "<<hit->stepLength()<<", "
143 <<", deposit energy: "<<hit->energyDeposit());
144 return hit;
145 }
147 const G4Step* hitStep) {
148 TrackHelper trkHelper{hitStep->GetTrack()};
151 if (m_writeHandle->empty() ||
152 m_writeHandle->back()->identify() != hitId ||
153 trkHelper.GenerateParticleLink() != m_writeHandle->back()->genParticleLink() ||
154 dec_G4TrkId(*m_writeHandle->back()) != hitStep->GetTrack()->GetTrackID()) {
155 return nullptr;
156 }
157 return m_writeHandle->back();
158 }
159}
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)
void setStore(AlignmentStorePtr store)
Adds the store to the Geometry context.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
const MuonGMR4::MuonDetectorManager * m_detMgr
Pointer to the underlying detector manager.
xAOD::MuonSimHit * lastSnapShot(const Identifier &gasGapId, const G4Step *hitStep)
Returns the last snap shot of the traversing particle.
SG::WriteHandle< xAOD::MuonSimHitContainer > m_writeHandle
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...
xAOD::MuonSimHit * propagateAndSaveStrip(const Identifier &hitId, const Amg::Transform3D &toGasGap, const G4Step *hitStep)
bool processStep(const G4Step *step) const
Checks whether the current step shall be processed at all.
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.
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
virtual bool isValid() override final
Can the handle be successfully dereferenced?
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.
const HepMcParticleLink & genParticleLink() const
Returns the link to the HepMC particle producing this hit.
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.
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