ATLAS Offline Software
Loading...
Searching...
No Matches
MuonG4R4::RpcSensitiveDetector Class Reference

Sensitive detector implementation to record G4 hits in the Rpc detectors. More...

#include <RpcSensitiveDetector.h>

Inheritance diagram for MuonG4R4::RpcSensitiveDetector:
Collaboration diagram for MuonG4R4::RpcSensitiveDetector:

Public Member Functions

 ~RpcSensitiveDetector ()=default
 Default destructor.
virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override final
 MuonSensitiveDetector (const std::string &name, const std::string &output_key, const std::string &trf_storeKey, const MuonGMR4::MuonDetectorManager *detMgr)
 Recycle the constructor from the MuonSensitiveDetector.
virtual void Initialize (G4HCofThisEvent *HCE) override final
 Create the output container at the beginning of the event.
bool msgLvl (const MSG::Level lvl) const
 Test the output level.
MsgStream & msg () const
 The standard message stream.
MsgStream & msg (const MSG::Level lvl) const
 The standard message stream.
void setLevel (MSG::Level lvl)
 Change the current logging level.

Protected Member Functions

bool processStep (const G4Step *step) const
 Checks whether the current step shall be processed at all.
ActsTrk::GeometryContext getGeoContext (const EventContext &ctx) const
 Returns the current geometry context in the event.
xAOD::MuonSimHitlastSnapShot (const Identifier &gasGapId, const G4Step *hitStep)
 Returns the last snap shot of the traversing particle.
xAOD::MuonSimHitpropagateAndSaveStrip (const Identifier &hitId, const Amg::Transform3D &toGasGap, const G4Step *hitStep)
 Records the G4Step in the sim hit.
xAOD::MuonSimHitsaveHit (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.

Protected Attributes

const MuonGMR4::MuonDetectorManagerm_detMgr {nullptr}
 Pointer to the underlying detector manager.

Private Member Functions

const MuonGMR4::RpcReadoutElementgetReadoutElement (const G4TouchableHistory *touchHist) const
 Retrieves the MuonReadoutElement associated to the rpc chamber in which the energy deposit is taking place.
Identifier getIdentifier (const ActsTrk::GeometryContext &gctx, const MuonGMR4::RpcReadoutElement *readOutEle, const Amg::Vector3D &hitAtGapPlane) const
 Identify the gas gap in which the G4 hit produced.
void initMessaging () const
 Initialize our message level and MessageSvc.

Private Attributes

std::string m_writeKey {}
 Key under which the output container is stored in the G4 event.
xAOD::MuonSimHitContainerm_outContainer {}
 Pointer to the MuonSimHit output container.
SG::ReadHandleKey< ActsTrk::DetectorAlignStorem_trfCacheKey
 ReadHandleKey to the DetectorAlignmentStore caching the relevant transformations needed in this event.
std::string m_nm
 Message source name.
boost::thread_specific_ptr< MsgStream > m_msg_tls
 MsgStream instance (a std::cout like with print-out levels).
std::atomic< IMessageSvc * > m_imsg { nullptr }
 MessageSvc pointer.
std::atomic< MSG::Level > m_lvl { MSG::NIL }
 Current logging level.
std::atomic_flag m_initialized ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
 Messaging initialized (initMessaging).

Detailed Description

Sensitive detector implementation to record G4 hits in the Rpc detectors.

The ProcessHits hook is called by Geant4 if the track enters a sensible Rpc gas gap volume. The TouchableHistory is used to deduce the associated readout element and then to identify the actual gas gap. The hit is then passed to the MuonSensitiveDetector class for event record

Definition at line 17 of file RpcSensitiveDetector.h.

Constructor & Destructor Documentation

◆ ~RpcSensitiveDetector()

MuonG4R4::RpcSensitiveDetector::~RpcSensitiveDetector ( )
default

Default destructor.

Member Function Documentation

◆ getGeoContext()

ActsTrk::GeometryContext MuonG4R4::MuonSensitiveDetector::getGeoContext ( const EventContext & ctx) const
protectedinherited

Returns the current geometry context in the event.

Parameters
ctxEvent context to access the store gate service

Definition at line 49 of file MuonSensitiveDetector.cxx.

49 {
50 ActsTrk::GeometryContext gctx{};
51 const ActsTrk::DetectorAlignStore* alignment{nullptr};
52 if (!SG::get(alignment, m_trfCacheKey, ctx).isSuccess()) {
53 THROW_EXCEPTION("Failed to retrieve "<<m_trfCacheKey.fullKey()<<".");
54 }
55 gctx.setStore(std::make_unique<DetectorAlignStore>(*alignment));
56 return gctx;
57 }
void setStore(AlignmentStorePtr store)
Adds the store to the Geometry context.
SG::ReadHandleKey< ActsTrk::DetectorAlignStore > m_trfCacheKey
ReadHandleKey to the DetectorAlignmentStore caching the relevant transformations needed in this event...
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10

◆ getIdentifier()

Identifier MuonG4R4::RpcSensitiveDetector::getIdentifier ( const ActsTrk::GeometryContext & gctx,
const MuonGMR4::RpcReadoutElement * readOutEle,
const Amg::Vector3D & hitAtGapPlane ) const
private

Identify the gas gap in which the G4 hit produced.

Parameters
gctxGeometry context to retrieve the center positions of the readout element's gas gaps
readOutEleThe previously identified readout element
hitAtGapPlanePosition of the G4 volume within the ATLAS coordinate system

Definition at line 66 of file RpcSensitiveDetector.cxx.

68 {
69 const RpcIdHelper& idHelper{m_detMgr->idHelperSvc()->rpcIdHelper()};
70 constexpr bool phiGap = false;
71 const Identifier firstChan = idHelper.channelID(readOutEle->identify(),
72 readOutEle->doubletZ(),
73 readOutEle->doubletPhi(), 1, phiGap, 1);
74
75 const Amg::Vector3D locHitPos{readOutEle->globalToLocalTransform(gctx, firstChan) *
76 hitAtGapPlane};
77 const double gapHalfWidth = readOutEle->stripEtaLength() / 2;
78 const double gapHalfLength = readOutEle->stripPhiLength()/ 2;
79 ATH_MSG_VERBOSE("Detector element: "<<m_detMgr->idHelperSvc()->toStringDetEl(firstChan)
80 <<" locPos: "<<Amg::toString(locHitPos, 2)
81 <<" gap thickness "<<readOutEle->gasGapPitch()
82 <<" gap width: "<<gapHalfWidth
83 <<" gap length: "<<gapHalfLength);
84 const int doubletPhi = std::abs(locHitPos.y()) > gapHalfWidth ? readOutEle->doubletPhiMax() :
85 readOutEle->doubletPhi();
86 const int gasGap = std::round(std::abs(locHitPos.z()) / readOutEle->gasGapPitch()) + 1;
87
88 return idHelper.channelID(readOutEle->identify(),
89 readOutEle->doubletZ(),
90 doubletPhi, gasGap, phiGap, 1);
91}
#define ATH_MSG_VERBOSE(x)
const MuonGMR4::MuonDetectorManager * m_detMgr
Pointer to the underlying detector manager.
Amg::Transform3D globalToLocalTransform(const ActsTrk::GeometryContext &ctx) const
Returns the transformation from the global ATLAS coordinate system into the local coordinate system o...
Identifier identify() const override final
Return the ATLAS identifier.
int doubletZ() const
Returns the doublet Z field of the MuonReadoutElement identifier.
int doubletPhi() const
Returns the doublet Phi field of the MuonReadoutElement identifier.
double stripPhiLength() const
Returns the length of a phi strip.
int doubletPhiMax() const
Returns the maximum phi panel.
double stripEtaLength() const
Returns the length of an eta strip.
double gasGapPitch() const
Returns the thickness of a RPC gasgap.
Identifier channelID(int stationName, int stationEta, int stationPhi, int doubletR, int doubletZ, int doubletPhi, int gasGap, int measuresPhi, int strip) const
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Eigen::Matrix< double, 3, 1 > Vector3D

◆ getReadoutElement()

const MuonGMR4::RpcReadoutElement * MuonG4R4::RpcSensitiveDetector::getReadoutElement ( const G4TouchableHistory * touchHist) const
private

Retrieves the MuonReadoutElement associated to the rpc chamber in which the energy deposit is taking place.

Parameters
touchHistTouchable history of the G4 track to find the proper volume name from which the readout element can be deduced

The fourth volume is the envelope volume of the rpc gas gap

Find the Detector element from the Identifier <STATIONETA>_<STATIONPHI>_<DOUBLETR>_<DOUBLETPHI>_<DOUBLETZ>

Definition at line 92 of file RpcSensitiveDetector.cxx.

92 {
94 const std::string stationVolume = touchHist->GetVolume(3)->GetName();
95
96 const std::vector<std::string> volumeTokens = tokenize(stationVolume, "_");
97 ATH_MSG_VERBOSE("Name of the station volume is "<<stationVolume);
98 if (volumeTokens.size() != 7) {
99 THROW_EXCEPTION(" Cannot deduce the station name from "<<stationVolume);
100 }
103 const std::string stName = volumeTokens[0].substr(0,3);
104 const int stationEta = atoi(volumeTokens[2]);
105 const int stationPhi = atoi(volumeTokens[3]);
106 const int doubletR = atoi(volumeTokens[4]);
107 const int doubletPhi = atoi(volumeTokens[5]);
108 const int doubletZ = atoi(volumeTokens[6]);
109 const RpcIdHelper& idHelper{m_detMgr->idHelperSvc()->rpcIdHelper()};
110
111 const Identifier detElId = idHelper.padID(idHelper.stationNameIndex(stName),
112 stationEta, stationPhi, doubletR, doubletZ, doubletPhi);
113 const RpcReadoutElement* readOutElem = m_detMgr->getRpcReadoutElement(detElId);
114 if (!readOutElem) {
115 THROW_EXCEPTION(" Failed to retrieve a valid detector element from "
116 <<m_detMgr->idHelperSvc()->toStringDetEl(detElId)<<" "<<stationVolume);
117 }
118 return readOutElem;
119}
int stationNameIndex(const std::string &name) const
Identifier padID(const Identifier &elementID, int doubletZ, int doubletPhi) const
std::vector< std::string > tokenize(std::string_view the_str, std::string_view delimiters)
Splits the string into smaller substrings.
int atoi(std::string_view str)
Helper functions to unpack numbers decoded in string into integers and doubles The strings are requir...
const std::string & stName(StIndex index)
convert StIndex into a string
constexpr uint8_t stationPhi
station Phi 1 to 8

◆ Initialize()

void MuonG4R4::MuonSensitiveDetector::Initialize ( G4HCofThisEvent * HCE)
finaloverridevirtualinherited

Create the output container at the beginning of the event.

Definition at line 38 of file MuonSensitiveDetector.cxx.

38 {
39 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
40 auto* hitVec = eventInfo->GetHitCollectionMap()->Find<MuonSimHitsVec>(m_writeKey);
41 if (!hitVec) {
42 THROW_EXCEPTION("The event does not contain a MuonSimHit container called "<<m_writeKey);
43 }
44 m_outContainer = hitVec->container.get();
45 } else {
46 THROW_EXCEPTION("There is no ATLAS event info");
47 }
48 }
static AtlasG4EventUserInfo * GetEventUserInfo()
std::string m_writeKey
Key under which the output container is stored in the G4 event.
xAOD::MuonSimHitContainer * m_outContainer
Pointer to the MuonSimHit output container.

◆ initMessaging()

void AthMessaging::initMessaging ( ) const
privateinherited

Initialize our message level and MessageSvc.

This method should only be called once.

Definition at line 39 of file AthMessaging.cxx.

40{
42 // If user did not set an explicit level, set a default
43 if (m_lvl == MSG::NIL) {
44 m_lvl = m_imsg ?
45 static_cast<MSG::Level>( m_imsg.load()->outputLevel(m_nm) ) :
46 MSG::INFO;
47 }
48}
std::string m_nm
Message source name.
std::atomic< IMessageSvc * > m_imsg
MessageSvc pointer.
std::atomic< MSG::Level > m_lvl
Current logging level.
IMessageSvc * getMessageSvc(bool quiet=false)

◆ lastSnapShot()

xAOD::MuonSimHit * MuonG4R4::MuonSensitiveDetector::lastSnapShot ( const Identifier & gasGapId,
const G4Step * hitStep )
protectedinherited

Returns the last snap shot of the traversing particle.

The G4 track must have stepped through the same volume. Otherwise, a nullptr is returned

Parameters
gasGapIdIdentifier of the gasGap to consider
hitStepPointer to the current step

There's only a snapshot if the last saved hit has the same Identifier & the same particle Link + G4Track Id

Definition at line 147 of file MuonSensitiveDetector.cxx.

148 {
149 TrackHelper trkHelper{hitStep->GetTrack()};
152 if (m_outContainer->empty() ||
153 m_outContainer->back()->identify() != hitId ||
154 trkHelper.GenerateParticleLink() != m_outContainer->back()->genParticleLink() ||
155 dec_G4TrkId(*m_outContainer->back()) != hitStep->GetTrack()->GetTrackID()) {
156 return nullptr;
157 }
158 return m_outContainer->back();
159 }
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52

◆ msg() [1/2]

MsgStream & AthMessaging::msg ( ) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 167 of file AthMessaging.h.

168{
169 MsgStream* ms = m_msg_tls.get();
170 if (!ms) {
171 if (!m_initialized.test_and_set()) initMessaging();
172 ms = new MsgStream(m_imsg,m_nm);
173 m_msg_tls.reset( ms );
174 }
175
176 ms->setLevel (m_lvl);
177 return *ms;
178}
boost::thread_specific_ptr< MsgStream > m_msg_tls
MsgStream instance (a std::cout like with print-out levels).
void initMessaging() const
Initialize our message level and MessageSvc.

◆ msg() [2/2]

MsgStream & AthMessaging::msg ( const MSG::Level lvl) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 182 of file AthMessaging.h.

183{ return msg() << lvl; }
MsgStream & msg() const
The standard message stream.

◆ msgLvl()

bool AthMessaging::msgLvl ( const MSG::Level lvl) const
inlineinherited

Test the output level.

Parameters
lvlThe message level to test against
Returns
boolean Indicating if messages at given level will be printed
Return values
trueMessages at level "lvl" will be printed

Definition at line 151 of file AthMessaging.h.

152{
153 // If user did not set explicit message level we have to initialize
154 // the messaging and retrieve the default via the MessageSvc.
155 if (m_lvl==MSG::NIL && !m_initialized.test_and_set()) initMessaging();
156
157 if (m_lvl <= lvl) {
158 msg() << lvl;
159 return true;
160 } else {
161 return false;
162 }
163}

◆ MuonSensitiveDetector()

MuonG4R4::MuonSensitiveDetector::MuonSensitiveDetector ( const std::string & name,
const std::string & output_key,
const std::string & trf_storeKey,
const MuonGMR4::MuonDetectorManager * detMgr )

Recycle the constructor from the MuonSensitiveDetector.

Definition at line 46 of file MuonSensitiveDetector.cxx.

30 :
31 G4VSensitiveDetector{name},
33 m_writeKey{output_key},
34 m_trfCacheKey{trfStore_key},
35 m_detMgr{detMgr} {
36 m_trfCacheKey.initialize().ignore();
37 }
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.

◆ ProcessHits()

G4bool MuonG4R4::RpcSensitiveDetector::ProcessHits ( G4Step * aStep,
G4TouchableHistory * ROhist )
finaloverridevirtual

The Rpc readout may have two strip panels per a single gasGap --> extrapolate the step towards the gap centre

Fetch the local -> global transformation

Definition at line 25 of file RpcSensitiveDetector.cxx.

25 {
26
27 if (!processStep(aStep)) {
28 return true;
29 }
30 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
31 const MuonGMR4::RpcReadoutElement* readOutEle = getReadoutElement(touchHist);
32 if (!readOutEle) {
33 return false;
34 }
35
37 if (!eventInfo) {
38 THROW_EXCEPTION("No AtlasG4EventUserInfo available");
39 }
40 const ActsTrk::GeometryContext gctx{getGeoContext(eventInfo->GetEventContext())};
41
42
43 const Amg::Transform3D localToGlobal = getTransform(touchHist, 0);
44 ATH_MSG_VERBOSE(" Track is inside volume "
45 <<touchHist->GetHistory()->GetTopVolume()->GetName()
46 <<" transformation: "<<Amg::toString(localToGlobal));
49 const Amg::Vector3D locPreStep{localToGlobal.inverse()*Amg::Hep3VectorToEigen(aStep->GetPreStepPoint()->GetPosition())};
50 const Amg::Vector3D locPostStep{localToGlobal.inverse()*Amg::Hep3VectorToEigen(aStep->GetPostStepPoint()->GetPosition())};
51
52 const Amg::Vector3D locStepDir = (locPostStep - locPreStep).unit();
53 const std::optional<double> lambda = Amg::intersect<3>(locPreStep, locStepDir, Amg::Vector3D::UnitX(), 0.);
54 Amg::Vector3D gapCentreCross = localToGlobal * ( (lambda ? 1. : 0.) *locPreStep + lambda.value_or(0.) * locStepDir);
55 const Identifier etaHitID = getIdentifier(gctx, readOutEle, gapCentreCross);
56 if (!etaHitID.is_valid()) {
57 ATH_MSG_VERBOSE("No valid hit found");
58 return true;
59 }
61 const Amg::Transform3D toGasGap{readOutEle->globalToLocalTransform(gctx, etaHitID)};
62 propagateAndSaveStrip(etaHitID, toGasGap, aStep);
63 return true;
64}
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
bool is_valid() const
Check if id is in a valid state.
ActsTrk::GeometryContext getGeoContext(const EventContext &ctx) const
Returns the current geometry context in the event.
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.
const MuonGMR4::RpcReadoutElement * getReadoutElement(const G4TouchableHistory *touchHist) const
Retrieves the MuonReadoutElement associated to the rpc chamber in which the energy deposit is taking ...
Identifier getIdentifier(const ActsTrk::GeometryContext &gctx, const MuonGMR4::RpcReadoutElement *readOutEle, const Amg::Vector3D &hitAtGapPlane) const
Identify the gas gap in which the G4 hit produced.
std::optional< double > intersect(const AmgVector(N)&posA, const AmgVector(N)&dirA, const AmgVector(N)&posB, const AmgVector(N)&dirB)
Calculates the point B' along the line B that's closest to a second line A.
Amg::Vector3D Hep3VectorToEigen(const CLHEP::Hep3Vector &CLHEPvector)
Converts a CLHEP-based CLHEP::Hep3Vector into an Eigen-based Amg::Vector3D.
Eigen::Affine3d Transform3D
Amg::Transform3D getTransform(const G4VTouchable *history, unsigned int level)
Extracts the local -> global transformation from a TouchableHistory at a given level.

◆ processStep()

bool MuonG4R4::MuonSensitiveDetector::processStep ( const G4Step * step) const
protectedinherited

Checks whether the current step shall be processed at all.

I.e. the particle needs to be charged, there's a minimum velocity needed and the step length must not vanish

Parameters
stepG4 step to consider

Reject secondary particles

Sensitive detector is only sensitive to charged particles or Geantinos

Definition at line 58 of file MuonSensitiveDetector.cxx.

58 {
59 const G4Track* currentTrack = aStep->GetTrack();
60 ATH_MSG_VERBOSE("Check whether step pdgId: "<<(*currentTrack)<<" will be processed.");
62 constexpr double velCutOff = 10.*Gaudi::Units::micrometer / Gaudi::Units::second;
63 if (aStep->GetStepLength() < std::numeric_limits<float>::epsilon() || currentTrack->GetVelocity() < velCutOff) {
64 ATH_MSG_VERBOSE("Step length is too short ");
65 return false;
66 }
68 if (currentTrack->GetDefinition()->GetPDGCharge() == 0.0) {
69 ATH_MSG_VERBOSE("Particle is neutral");
70 return currentTrack->GetDefinition() == G4Geantino::GeantinoDefinition() ||
71 currentTrack->GetDefinition() == G4ChargedGeantino::ChargedGeantinoDefinition();
72 }
73 return true;
74 }

◆ propagateAndSaveStrip()

xAOD::MuonSimHit * MuonG4R4::MuonSensitiveDetector::propagateAndSaveStrip ( const Identifier & hitId,
const Amg::Transform3D & toGasGap,
const G4Step * hitStep )
protectedinherited

Records the G4Step in the sim hit.

Hits are usually expressed at the center point of the step. Except for low-energy electrons where the endpoint of the random walk path is just extended

Parameters
hitIdIdentifier of the gas gap or the tube in which the hit is recorded
toGaGapTransform from the ATLAS global coordinates into the sensitive volume
hitStepStep from which the particle's information and the pre and post step positions are fetched

Fetch the step end-points

Hit direction as the momentum direction

Electrons randomly work through the gas instead of drifting through the gas -> take the prestep of the first snap shop inside the gas as pre step point

Definition at line 75 of file MuonSensitiveDetector.cxx.

77 {
78
79 const G4Track* currentTrack = aStep->GetTrack();
81 const Amg::Vector3D locPostStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPostStepPoint()->GetPosition())};
82 Amg::Vector3D locPreStep{toGasGap*Amg::Hep3VectorToEigen(aStep->GetPreStepPoint()->GetPosition())};
83
85 Amg::Vector3D locHitDir = toGasGap.linear() * Amg::Hep3VectorToEigen(currentTrack->GetMomentumDirection());
86
89 xAOD::MuonSimHit* prevHit = lastSnapShot(hitID, aStep);
90 if (std::abs(currentTrack->GetParticleDefinition()->GetPDGEncoding()) == 11 && prevHit) {
91 locPreStep = xAOD::toEigen(prevHit->localPosition()) - 0.5* prevHit->stepLength() *
92 xAOD::toEigen(prevHit->localDirection());
93
94 locHitDir = (locPostStep - locPreStep).unit();
95 }
96 const Amg::Vector3D locHitPos = 0.5* (locPreStep + locPostStep);
97 ATH_MSG_VERBOSE( m_detMgr->idHelperSvc()->toStringGasGap(hitID)<<" - track: "<<(*currentTrack)
98 <<", deposit: "<<aStep->GetTotalEnergyDeposit()<<", -- local coords: "
99 <<"prestep: "<<Amg::toString(locPreStep)<<", post step: "<<Amg::toString(locPostStep)
100 <<" mid point: "<< Amg::toString(locHitPos)<<", direction: "<<Amg::toString(locHitDir)
101 <<", deposit: "<<aStep->GetTotalEnergyDeposit());
102
103 const double globalTime = currentTrack->GetGlobalTime() + locHitDir.dot(locPostStep - locHitPos) / currentTrack->GetVelocity();
104
105 xAOD::MuonSimHit* newHit = saveHit(hitID, locHitPos, locHitDir, globalTime, aStep);
106 if (prevHit) {
107 newHit->setStepLength((locPostStep - locPreStep).mag());
108 }
109 return newHit;
110 }
Scalar mag() const
mag method
xAOD::MuonSimHit * lastSnapShot(const Identifier &gasGapId, const G4Step *hitStep)
Returns the last snap shot of the traversing particle.
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.
ConstVectorMap< 3 > localDirection() const
Returns the local direction of the traversing particle.
float stepLength() const
Returns the path length of the G4 step.
void setStepLength(const float length)
Set the path length of the G4 step.
ConstVectorMap< 3 > localPosition() const
Returns the local postion of the traversing particle.
MuonSimHit_v1 MuonSimHit
Defined the version of the MuonSimHit.
Definition MuonSimHit.h:12

◆ saveHit()

xAOD::MuonSimHit * MuonG4R4::MuonSensitiveDetector::saveHit ( const Identifier & hitId,
const Amg::Vector3D & hitPos,
const Amg::Vector3D & hitDir,
const double globTime,
const G4Step * hitStep )
protectedinherited

Saves the current Step as a xAOD::MuonSimHit snapshot.

Parameters
hitIdIdentifier of the gasGap/ tube where the hit was deposited
hitPosLocal position of the hit expressed w.r.t gas gap coordinate system
hitDirLocal direction of the hit expressed w.r.t gas gap coordinate system
globTimeGlobal time of the hit
hitStepPointer to the step from which the hit information was derived.

Definition at line 111 of file MuonSensitiveDetector.cxx.

115 {
116 const G4Track* currentTrack = aStep->GetTrack();
117 TrackHelper trHelper{currentTrack};
118 // If Geant4 propagates the same track through the same volume just update the last hit and don't write a new one
119 xAOD::MuonSimHit* hit = lastSnapShot(hitId, aStep);
120 bool newHit{false};
121 if (!hit) {
122 hit = m_outContainer->push_back(std::make_unique<xAOD::MuonSimHit>());
123 newHit = true;
124 }
125 dec_G4TrkId(*hit) = currentTrack->GetTrackID();
126 hit->setIdentifier(hitId);
127 hit->setLocalPosition(xAOD::toStorage(hitPos));
129 hit->setMass(currentTrack->GetDefinition()->GetPDGMass());
130 hit->setGlobalTime(globTime);
131 hit->setPdgId(currentTrack->GetDefinition()->GetPDGEncoding());
132 hit->setEnergyDeposit(aStep->GetTotalEnergyDeposit() + (newHit ? 0. : hit->energyDeposit()));
133 hit->setKineticEnergy(currentTrack->GetKineticEnergy());
135 hit->setStepLength(aStep->GetStepLength());
136
137 ATH_MSG_VERBOSE("Save new hit "<<m_detMgr->idHelperSvc()->toString(hitId)
138 <<", pdgId: "<<hit->pdgId()
139 <<", "<<hit->genParticleLink()
140 <<", trackId: "<<currentTrack->GetTrackID()<<", "
141 <<", "<<hit->genParticleLink().cptr()<<std::endl
142 <<"pos: "<<Amg::toString(hitPos)<<", dir: "<<Amg::toString(hitDir)<<", time: "<<globTime
143 <<", energy: "<<hit->kineticEnergy()<<", stepLength: "<<hit->stepLength()<<", "
144 <<", deposit energy: "<<hit->energyDeposit());
145 return hit;
146 }
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
int pdgId() const
Returns the pdgID of the traversing particle.
void setIdentifier(const Identifier &id)
Sets the global ATLAS identifier.
float energyDeposit() const
Returns the energy deposited by the traversing particle inside the gas volume.
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.
MeasVector< N > toStorage(const AmgVector(N)&amgVec)
Converts the double precision of the AmgVector into the floating point storage precision of the MeasV...

◆ setLevel()

void AthMessaging::setLevel ( MSG::Level lvl)
inherited

Change the current logging level.

Use this rather than msg().setLevel() for proper operation with MT.

Definition at line 28 of file AthMessaging.cxx.

29{
30 m_lvl = lvl;
31}

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
mutableprivateinherited

Messaging initialized (initMessaging).

Definition at line 141 of file AthMessaging.h.

◆ m_detMgr

const MuonGMR4::MuonDetectorManager* MuonG4R4::MuonSensitiveDetector::m_detMgr {nullptr}
protectedinherited

Pointer to the underlying detector manager.

Definition at line 107 of file MuonSensitiveDetector.h.

107{nullptr};

◆ m_imsg

std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr }
mutableprivateinherited

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_lvl

std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL }
mutableprivateinherited

Current logging level.

Definition at line 138 of file AthMessaging.h.

138{ MSG::NIL };

◆ m_msg_tls

boost::thread_specific_ptr<MsgStream> AthMessaging::m_msg_tls
mutableprivateinherited

MsgStream instance (a std::cout like with print-out levels).

Definition at line 132 of file AthMessaging.h.

◆ m_nm

std::string AthMessaging::m_nm
privateinherited

Message source name.

Definition at line 129 of file AthMessaging.h.

◆ m_outContainer

xAOD::MuonSimHitContainer* MuonG4R4::MuonSensitiveDetector::m_outContainer {}
privateinherited

Pointer to the MuonSimHit output container.

Definition at line 60 of file MuonSensitiveDetector.h.

60{};

◆ m_trfCacheKey

SG::ReadHandleKey<ActsTrk::DetectorAlignStore> MuonG4R4::MuonSensitiveDetector::m_trfCacheKey
privateinherited

ReadHandleKey to the DetectorAlignmentStore caching the relevant transformations needed in this event.

Definition at line 63 of file MuonSensitiveDetector.h.

◆ m_writeKey

std::string MuonG4R4::MuonSensitiveDetector::m_writeKey {}
privateinherited

Key under which the output container is stored in the G4 event.

Definition at line 58 of file MuonSensitiveDetector.h.

58{};

The documentation for this class was generated from the following files: