ATLAS Offline Software
Loading...
Searching...
No Matches
RpcSensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
7#include "G4ThreeVector.hh"
8#include "G4Trd.hh"
9#include "G4Geantino.hh"
10#include "G4ChargedGeantino.hh"
11
12#include "MCTruth/TrackHelper.h"
13#include <sstream>
14
17#include "GaudiKernel/SystemOfUnits.h"
18#include "GeoModelKernel/throwExcept.h"
19
20using namespace MuonGMR4;
21using namespace CxxUtils;
22using namespace ActsTrk;
23
24// construction/destruction
25namespace MuonG4R4 {
26
27
28G4bool RpcSensitiveDetector::ProcessHits(G4Step* aStep,G4TouchableHistory*) {
29
30 if (!processStep(aStep)) {
31 return true;
32 }
33 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
34 const MuonGMR4::RpcReadoutElement* readOutEle = getReadoutElement(touchHist);
35 if (!readOutEle) {
36 return false;
37 }
38
40
41 const Amg::Transform3D localToGlobal = getTransform(touchHist, 0);
42 ATH_MSG_VERBOSE(" Track is inside volume "
43 <<touchHist->GetHistory()->GetTopVolume()->GetName()
44 <<" transformation: "<<Amg::toString(localToGlobal));
47 const Amg::Vector3D locPreStep{localToGlobal.inverse()*Amg::Hep3VectorToEigen(aStep->GetPreStepPoint()->GetPosition())};
48 const Amg::Vector3D locPostStep{localToGlobal.inverse()*Amg::Hep3VectorToEigen(aStep->GetPostStepPoint()->GetPosition())};
49
50 const Amg::Vector3D locStepDir = (locPostStep - locPreStep).unit();
51 const std::optional<double> lambda = Amg::intersect<3>(locPreStep, locStepDir, Amg::Vector3D::UnitX(), 0.);
52 Amg::Vector3D gapCentreCross = localToGlobal * ( (lambda ? 1. : 0.) *locPreStep + lambda.value_or(0.) * locStepDir);
53 const Identifier etaHitID = getIdentifier(gctx, readOutEle, gapCentreCross, false);
54 if (!etaHitID.is_valid()) {
55 ATH_MSG_VERBOSE("No valid hit found");
56 return true;
57 }
59 const Amg::Transform3D toGasGap{readOutEle->globalToLocalTrans(gctx, etaHitID)};
60 propagateAndSaveStrip(etaHitID, toGasGap, aStep);
61 return true;
62}
63
65 const MuonGMR4::RpcReadoutElement* readOutEle,
66 const Amg::Vector3D& hitAtGapPlane, bool phiGap) const {
67 const RpcIdHelper& idHelper{m_detMgr->idHelperSvc()->rpcIdHelper()};
68
69 const Identifier firstChan = idHelper.channelID(readOutEle->identify(),
70 readOutEle->doubletZ(),
71 readOutEle->doubletPhi(), 1, phiGap, 1);
72
73 const Amg::Vector3D locHitPos{readOutEle->globalToLocalTrans(gctx, firstChan) *
74 hitAtGapPlane};
75 const double gapHalfWidth = readOutEle->stripEtaLength() / 2;
76 const double gapHalfLength = readOutEle->stripPhiLength()/ 2;
77 ATH_MSG_VERBOSE("Detector element: "<<m_detMgr->idHelperSvc()->toStringDetEl(firstChan)
78 <<" locPos: "<<Amg::toString(locHitPos, 2)
79 <<" gap thickness "<<readOutEle->gasGapPitch()
80 <<" gap width: "<<gapHalfWidth
81 <<" gap length: "<<gapHalfLength);
82 const int doubletPhi = std::abs(locHitPos.y()) > gapHalfWidth ? readOutEle->doubletPhiMax() :
83 readOutEle->doubletPhi();
84 const int gasGap = std::round(std::abs(locHitPos.z()) / readOutEle->gasGapPitch()) + 1;
85
86 return idHelper.channelID(readOutEle->identify(),
87 readOutEle->doubletZ(),
88 doubletPhi, gasGap, phiGap, 1);
89}
90const MuonGMR4::RpcReadoutElement* RpcSensitiveDetector::getReadoutElement(const G4TouchableHistory* touchHist) const {
92 const std::string stationVolume = touchHist->GetVolume(3)->GetName();
93
94 const std::vector<std::string> volumeTokens = tokenize(stationVolume, "_");
95 ATH_MSG_VERBOSE("Name of the station volume is "<<stationVolume);
96 if (volumeTokens.size() != 7) {
97 THROW_EXCEPTION(" Cannot deduce the station name from "<<stationVolume);
98 }
101 const std::string stName = volumeTokens[0].substr(0,3);
102 const int stationEta = atoi(volumeTokens[2]);
103 const int stationPhi = atoi(volumeTokens[3]);
104 const int doubletR = atoi(volumeTokens[4]);
105 const int doubletPhi = atoi(volumeTokens[5]);
106 const int doubletZ = atoi(volumeTokens[6]);
107 const RpcIdHelper& idHelper{m_detMgr->idHelperSvc()->rpcIdHelper()};
108
109 const Identifier detElId = idHelper.padID(idHelper.stationNameIndex(stName),
110 stationEta, stationPhi, doubletR, doubletZ, doubletPhi);
111 const RpcReadoutElement* readOutElem = m_detMgr->getRpcReadoutElement(detElId);
112 if (!readOutElem) {
113 THROW_EXCEPTION(" Failed to retrieve a valid detector element from "
114 <<m_detMgr->idHelperSvc()->toStringDetEl(detElId)<<" "<<stationVolume);
115 }
116 return readOutElem;
117}
118}
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
#define ATH_MSG_VERBOSE(x)
bool is_valid() const
Check if id is in a valid state.
const MuonGMR4::MuonDetectorManager * m_detMgr
Pointer to the underlying detector manager.
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.
ActsTrk::GeometryContext getGeoContext() const
Returns the current geometry context in the event.
const MuonGMR4::RpcReadoutElement * getReadoutElement(const G4TouchableHistory *touchHist) const
Retrieves the matching readout element to a G4 hit.
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override final
Identifier getIdentifier(const ActsTrk::GeometryContext &gctx, const MuonGMR4::RpcReadoutElement *readOutEle, const Amg::Vector3D &hitAtGapPlane, bool phiGap) const
Extracts the gasGap Identifier of the hit expressed at the origin of the local gasGap system.
Identifier identify() const override final
Return the athena identifier.
Amg::Transform3D globalToLocalTrans(const ActsTrk::GeometryContext &ctx) const
Transformations to translate between local <-> global coordinates.
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.
int stationNameIndex(const std::string &name) const
Identifier padID(const Identifier &elementID, int doubletZ, int doubletPhi) const
Identifier channelID(int stationName, int stationEta, int stationPhi, int doubletR, int doubletZ, int doubletPhi, int gasGap, int measuresPhi, int strip) const
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
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.
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
std::vector< std::string > tokenize(const std::string &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...
Include the common definitions from the MuonReadoutGeometry.
Amg::Transform3D getTransform(const G4VTouchable *history, unsigned int level)
Extracts the local -> global transformation from a TouchableHistory at a given level.
The ReadoutGeomCnvAlg converts the Run4 Readout geometry build from the GeoModelXML into the legacy M...
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10