ATLAS Offline Software
Loading...
Searching...
No Matches
MuonDetDescr/MuonGeoModelTest/src/GeoModelRpcTest.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4#include "GeoModelRpcTest.h"
5
6#include <fstream>
7#include <iostream>
8
13#include "GaudiKernel/SystemOfUnits.h"
14
15
16namespace MuonGM {
17
19 ATH_CHECK(m_tree.write());
20 return StatusCode::SUCCESS;
21}
23 ATH_CHECK(m_detMgrKey.initialize());
24 ATH_CHECK(m_idHelperSvc.retrieve());
25 ATH_CHECK(m_tree.init(this));
26 const RpcIdHelper& idHelper{m_idHelperSvc->rpcIdHelper()};
27 auto translateTokenList = [this, &idHelper](const std::vector<std::string>& chNames){
28
29 std::set<Identifier> transcriptedIds{};
30 for (const std::string& token : chNames) {
31 if (token.size() != 6) {
32 ATH_MSG_WARNING("Wrong format given for "<<token<<". Expecting 6 characters");
33 continue;
34 }
36 const std::string statName = token.substr(0, 3);
37 const unsigned statEta = std::atoi(token.substr(3, 1).c_str()) * (token[4] == 'A' ? 1 : -1);
38 const unsigned statPhi = std::atoi(token.substr(5, 1).c_str());
39 bool isValid{false};
40 const Identifier eleId = idHelper.elementID(statName, statEta, statPhi, 1, isValid);
41 if (!isValid) {
42 ATH_MSG_WARNING("Failed to deduce a station name for " << token);
43 continue;
44 }
45 transcriptedIds.insert(eleId);
46 std::copy_if(idHelper.detectorElement_begin(), idHelper.detectorElement_end(),
47 std::inserter(transcriptedIds, transcriptedIds.end()),
48 [&eleId, &idHelper](const Identifier& copyMe){
49 return idHelper.stationName(copyMe) == idHelper.stationName(eleId) &&
50 idHelper.stationEta(copyMe) == idHelper.stationEta(eleId) &&
51 idHelper.stationPhi(copyMe) == idHelper.stationPhi(eleId);
52 });
53 }
54 return transcriptedIds;
55 };
56
57 std::vector <std::string>& selectedSt = m_selectStat.value();
58 const std::vector <std::string>& excludedSt = m_excludeStat.value();
59 selectedSt.erase(std::remove_if(selectedSt.begin(), selectedSt.end(),
60 [&excludedSt](const std::string& token){
61 return std::ranges::find(excludedSt, token) != excludedSt.end();
62 }), selectedSt.end());
63
64 if (selectedSt.size()) {
65 m_testStations = translateTokenList(selectedSt);
66 std::stringstream sstr{};
67 for (const Identifier& id : m_testStations) {
68 sstr<<" *** "<<m_idHelperSvc->toString(id)<<std::endl;
69 }
70 ATH_MSG_INFO("Test only the following stations "<<std::endl<<sstr.str());
71 } else {
72 const std::set<Identifier> excluded = translateTokenList(excludedSt);
74 for(auto itr = idHelper.detectorElement_begin();
75 itr!= idHelper.detectorElement_end();++itr){
76 if (!excluded.count(*itr)) {
77 m_testStations.insert(*itr);
78 }
79 }
81 if (!excluded.empty()) {
82 std::stringstream excluded_report{};
83 for (const Identifier& id : excluded){
84 excluded_report << " *** " << m_idHelperSvc->toStringDetEl(id) << std::endl;
85 }
86 ATH_MSG_INFO("Test all station except the following excluded ones " << std::endl << excluded_report.str());
87 }
88 }
89 return StatusCode::SUCCESS;
90}
92 const EventContext& ctx{Gaudi::Hive::currentContext()};
94 if (!detMgr.isValid()) {
95 ATH_MSG_FATAL("Failed to retrieve MuonDetectorManager "
96 << m_detMgrKey.fullKey());
97 return StatusCode::FAILURE;
98 }
99 for (const Identifier& test_me : m_testStations) {
100 ATH_MSG_VERBOSE("Test retrieval of Mdt detector element "
101 << m_idHelperSvc->toStringDetEl(test_me));
102 const RpcReadoutElement* reElement = detMgr->getRpcReadoutElement(test_me);
103 if (!reElement) {
104 ATH_MSG_VERBOSE("Detector element is invalid");
105 continue;
106 }
108 if (reElement->identify() != test_me) {
109 ATH_MSG_FATAL("Expected to retrieve "
110 << m_idHelperSvc->toStringDetEl(test_me) << ". But got instead "
111 << m_idHelperSvc->toStringDetEl(reElement->identify()));
112 return StatusCode::FAILURE;
113 }
114 ATH_CHECK(dumpToTree(ctx, reElement));
115 if (m_idHelperSvc->stationNameString(reElement->identify()) == "BIS") continue;
116 const RpcIdHelper& idHelper{m_idHelperSvc->rpcIdHelper()};
117 for (int gasGap = 1 ; gasGap <= reElement->numberOfLayers(); ++gasGap) {
118 for (int doubPhi = reElement->getDoubletPhi(); doubPhi <= reElement->NphiStripPanels(); ++doubPhi) {
119 for (bool measPhi: {false, true}) {
120 for (int strip = 1 ; strip < reElement->Nstrips(measPhi); ++strip) {
121 bool isValid{false};
122 const Identifier stripId = idHelper.channelID(test_me,
123 idHelper.doubletZ(test_me),
124 doubPhi,
125 gasGap, measPhi, strip, isValid);
126 if (!isValid) {
127 ATH_MSG_DEBUG("Could not construct Identifier from "
128 <<m_idHelperSvc->toStringChamber(test_me)
129 <<", gasGap:"<<gasGap<<", doubletPhi: "<<doubPhi
130 <<", measurePhi "<<measPhi<<", strip: "<<strip);
131 continue;
132 }
133 Amg::Vector2D locPos2D{Amg::Vector2D::Zero()};
134 if (!reElement->stripPosition(stripId, locPos2D)) {
135 ATH_MSG_FATAL("Could not retrieve the local strip position for "<<m_idHelperSvc->toString(stripId));
136 return StatusCode::FAILURE;
137 }
138 const Trk::Surface& planeSurf{reElement->surface(stripId)};
139 const Amg::Vector3D globPos3D = planeSurf.localToGlobal(locPos2D);
141 const Amg::Vector3D stripPos = reElement->stripPos(stripId);
142 if ( (stripPos - globPos3D).mag() > std::numeric_limits<float>::epsilon()) {
143 ATH_MSG_FATAL("Retrieving the strip position in two different ways leads to two distinct results "
144 <<Amg::toString(stripPos)<<" as reference. Second path "<<Amg::toString(globPos3D));
145 return StatusCode::FAILURE;
146 }
147 const Trk::Surface& laySurf{reElement->surface(stripId)};
148 const double stripLen = 0.5 *reElement->StripLength(measPhi) - 1. * Gaudi::Units::cm;
149 if (!laySurf.insideBounds(locPos2D)){
150 ATH_MSG_FATAL("The strip center "<<Amg::toString(locPos2D)<<" of "<<m_idHelperSvc->toString(stripId)
151 <<" is outside bounds "<<laySurf.bounds());
152 return StatusCode::FAILURE;
153 }
154 if (!laySurf.insideBounds(locPos2D + stripLen * Amg::Vector2D::UnitY())){
155 ATH_MSG_FATAL("The right strip edge "<<Amg::toString(locPos2D + stripLen * Amg::Vector2D::UnitY())<<
156 " of "<<m_idHelperSvc->toString(stripId)<<" is outside bounds "<<laySurf.bounds());
157 return StatusCode::FAILURE;
158 }if (!laySurf.insideBounds(locPos2D - stripLen * Amg::Vector2D::UnitY())){
159 ATH_MSG_FATAL("The left strip edge "<<Amg::toString(locPos2D - stripLen * Amg::Vector2D::UnitY())
160 <<" of "<<m_idHelperSvc->toString(stripId)<<" is outside bounds "<<laySurf.bounds());
161 return StatusCode::FAILURE;
162 }
163 }
164 }
165 }
166 }
167 }
168 return StatusCode::SUCCESS;
169}
170StatusCode GeoModelRpcTest::dumpToTree(const EventContext& ctx, const RpcReadoutElement* readoutEle) {
171
172 m_stIndex = readoutEle->getStationIndex();
173 m_stEta = readoutEle->getStationEta();
174 m_stPhi = readoutEle->getStationPhi();
175 m_doubletR = readoutEle->getDoubletR();
176 m_doubletZ = readoutEle->getDoubletZ();
177 m_doubletPhi = readoutEle->getDoubletPhi();
178 m_chamberDesign = readoutEle->getTechnologyName();
179
180 m_numStripsEta = readoutEle->Nstrips(false);
181 m_numStripsPhi = readoutEle->Nstrips(true);
182 m_numRpcLayers = readoutEle->numberOfLayers();
183 m_numPhiPanels = readoutEle->NphiStripPanels();
184
185 m_stripEtaPitch = readoutEle->StripPitch(false);
186 m_stripPhiPitch = readoutEle->StripPitch(true);
187 m_stripEtaWidth = readoutEle->StripWidth(false);
188 m_stripPhiWidth = readoutEle->StripWidth(true);
189 m_stripEtaLength = readoutEle->StripLength(false);
190 m_stripPhiLength = readoutEle->StripLength(true);
191
192
193 const RpcIdHelper& idHelper{m_idHelperSvc->rpcIdHelper()};
194
195 const Amg::Transform3D& trans{readoutEle->absTransform()};
196 m_readoutTransform = trans;
197 const MuonGM::MuonStation* station = readoutEle->parentMuonStation();
198 m_alignableNode = station->getGeoTransform()->getDefTransform() *
199 station->getNativeToAmdbLRS().inverse();
200 if (station->hasALines()) {
201 m_ALineTransS = station->getALine_tras();
202 m_ALineTransT = station->getALine_traz();
203 m_ALineTransZ = station->getALine_trat();
204 m_ALineRotS = station->getALine_rots();
205 m_ALineRotT = station->getALine_rotz();
206 m_ALineRotZ = station->getALine_rott();
207 }
208 const int maxDoubPhi = std::max(readoutEle->getDoubletPhi(), readoutEle->NphiStripPanels());
209 for (int doubPhi = readoutEle->getDoubletPhi(); doubPhi <= maxDoubPhi; ++doubPhi) {
210 for (int gap = 1; gap <= readoutEle->numberOfLayers(); ++gap) {
211 for (bool measPhi : {false, true}) {
212 unsigned int numStrip = readoutEle->Nstrips(measPhi);
213 for (unsigned int strip = 1; strip <= numStrip ; ++strip) {
214 bool isValid{false};
215 const Identifier stripID = idHelper.channelID(readoutEle->identify(),
216 readoutEle->getDoubletZ(),
217 doubPhi, gap, measPhi, strip, isValid);
218 if (!isValid) {
219 ATH_MSG_WARNING("Invalid Identifier detected for readout element "
220 <<m_idHelperSvc->toStringDetEl(readoutEle->identify())
221 <<" gap: "<<gap<<" strip: "<<strip<<" meas phi: "<<measPhi);
222 continue;
223 }
224 Amg::Vector2D lStripPos{Amg::Vector2D::Zero()};
225 if (!readoutEle->stripPosition(stripID, lStripPos)){
226 ATH_MSG_FATAL("Failed to obtain strip position for "<<m_idHelperSvc->toString(stripID));
227 return StatusCode::FAILURE;
228 }
229
230 m_stripPos.push_back(readoutEle->stripPos(stripID));
231 m_locPos.push_back(lStripPos);
232 m_stripPosGasGap.push_back(gap);
233 m_stripPosMeasPhi.push_back(measPhi);
234 m_stripPosNum.push_back(strip);
235 m_stripDblPhi.push_back(doubPhi);
236
237 if (strip != 1) continue;
238 const Amg::Transform3D locToGlob = readoutEle->transform(stripID);
239 m_stripRot.push_back(locToGlob);
240 m_stripRotGasGap.push_back(gap);
241 m_stripRotMeasPhi.push_back(measPhi);
242 m_stripRotDblPhi.push_back(doubPhi);
243
244 }
245 }
246 }
247 }
248
249 return m_tree.fill(ctx) ? StatusCode::SUCCESS : StatusCode::FAILURE;
250}
251
252
253}
Scalar mag() const
mag method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
MuonVal::ScalarBranch< float > & m_stripEtaPitch
Strip dimensions.
std::set< Identifier > m_testStations
Set of stations to be tested.
MuonVal::ThreeVectorBranch m_stripPos
Strip positions.
MuonVal::ScalarBranch< uint8_t > & m_numRpcLayers
Number of eta & phi gas gaps.
MuonVal::ScalarBranch< float > & m_ALineTransS
Alignment parameters.
Gaudi::Property< std::vector< std::string > > m_selectStat
String should be formated like <stationName><stationEta><A/C><stationPhi>
Gaudi::Property< std::vector< std::string > > m_excludeStat
MuonVal::ScalarBranch< std::string > & m_chamberDesign
MuonVal::ScalarBranch< unsigned short > & m_stIndex
Identifier of the readout element.
MuonVal::ScalarBranch< uint8_t > & m_numStripsEta
Number of strips, strip pitch in eta & phi direction.
SG::ReadCondHandleKey< MuonGM::MuonDetectorManager > m_detMgrKey
MuonDetectorManager from the conditions store.
MuonVal::MuonTesterTree m_tree
Write a TTree for validation purposes.
MuonVal::CoordTransformBranch m_readoutTransform
Transformation of the readout element (Translation, ColX, ColY, ColZ)
MuonVal::CoordSystemsBranch m_stripRot
Rotation matrix of the respective layers.
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
StatusCode dumpToTree(const EventContext &ctx, const RpcReadoutElement *readoutEle)
virtual const Trk::PlaneSurface & surface() const override
access to chamber surface (phi orientation), uses the first gas gap
virtual const Amg::Transform3D & transform() const override
Return local to global transform.
Identifier identify() const override final
Returns the ATLAS Identifier of the MuonReadOutElement.
double getALine_trat() const
const Amg::Transform3D & getNativeToAmdbLRS() const
double getALine_rott() const
double getALine_rots() const
const GeoAlignableTransform * getGeoTransform() const
bool hasALines() const
double getALine_tras() const
double getALine_rotz() const
double getALine_traz() const
An RpcReadoutElement corresponds to a single RPC module; therefore typicaly a barrel muon station con...
int getDoubletPhi() const
return DoubletPhi value for the given readout element, be aware that one RE can contain two DoubletPh...
virtual bool stripPosition(const Identifier &id, Amg::Vector2D &pos) const override final
strip position If the strip number is outside the range of valid strips, the function will return fal...
int NphiStripPanels() const
returns the number of phi strip panels (1 or 2)
int getDoubletR() const
return DoubletR value for the given readout element
double StripLength(bool measphi) const
returns the strip length for the phi or eta plane
int getDoubletZ() const
return DoubletZ value for the given readout element
int Nstrips(bool measphi) const
returns the number of strips for the phi or eta plane
virtual int numberOfLayers(bool measphi=true) const override final
number of layers in phi/eta projection, same for eta/phi planes
double StripWidth(bool measphi) const
returns the strip width for the phi or eta plane
double StripPitch(bool measphi) const
returns the strip pitch for the phi or eta plane
const_id_iterator detectorElement_begin() const
Iterators over full set of ids.
const_id_iterator detectorElement_end() const
Identifier channelID(int stationName, int stationEta, int stationPhi, int doubletR, int doubletZ, int doubletPhi, int gasGap, int measuresPhi, int strip) const
Identifier elementID(int stationName, int stationEta, int stationPhi, int doubletR) const
int doubletZ(const Identifier &id) const
Abstract Base Class for tracking surfaces.
virtual bool insideBounds(const Amg::Vector2D &locpos, double tol1=0., double tol2=0.) const =0
virtual methods to be overwritten by the inherited surfaces
virtual const SurfaceBounds & bounds() const =0
Surface Bounds method.
virtual void localToGlobal(const Amg::Vector2D &locp, const Amg::Vector3D &mom, Amg::Vector3D &glob) const =0
Specified by each surface type: LocalToGlobal method without dynamic memory allocation.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the Athena extensions are properly loaded.
Definition GeoMuonHits.h:27
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.