ATLAS Offline Software
Loading...
Searching...
No Matches
RPCDQUtils.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// C/C++
6#include <math.h>
7
8#include <fstream>
9
10// Athena
15// Local
16
17#include "RPCDQUtils.h"
18
19//========================================================================================================
20// struct RpcPanel
21//========================================================================================================
22RpcPanel::RpcPanel(Identifier id, const RpcIdHelper &rpcIdHelper) {
23 FillRpcId(id, rpcIdHelper);
24}
25
26//========================================================================================================
27void RpcPanel::FillRpcId(Identifier id, const RpcIdHelper &rpcIdHelper) {
28 //
29 // Read RPC data
30 //
31 stationName = rpcIdHelper.stationName(id);
32 stationEta = rpcIdHelper.stationEta(id);
33 stationPhi = rpcIdHelper.stationPhi(id);
34 doubletR = rpcIdHelper.doubletR(id);
35 doubletZ = rpcIdHelper.doubletZ(id);
36 doubletPhi = rpcIdHelper.doubletPhi(id);
37 gasGap = rpcIdHelper.gasGap(id);
38 measPhi = rpcIdHelper.measuresPhi(id);
39
40 panelId =
43
44 // panel_index = stationName*6144 + stationPhi*768 + (stationEta+9)*48 +
45 // doubletR*24 + doubletPhi*12 + doubletZ*4 + gasGap*2 + measPhi;
46}
47
48//========================================================================================================
50 const MuonGM::RpcReadoutElement *readoutEl_,
51 const int doubletZ_, const int doubletPhi_,
52 const int gasgap_, const int measPhi_)
53 : readoutEl(readoutEl_),
54 doubletR(readoutEl_->getDoubletR()),
55 doubletZ(doubletZ_),
56 doubletPhi(doubletPhi_),
57 gasGap(gasgap_),
58 measPhi(measPhi_) {
59 const RpcIdHelper &rpcIdHelper = idHelperSvc.rpcIdHelper();
60 const Identifier readEl_id = readoutEl->identify();
61
62 stationName = rpcIdHelper.stationName(readEl_id);
64
65 stationEta = rpcIdHelper.stationEta(readEl_id);
66 stationPhi = rpcIdHelper.stationPhi(readEl_id);
67
68 //
69 // Get Identifier for this panel
70 //
71 panelId =
74 if (panel_valid) {
75 panel_str = idHelperSvc.toString(panelId);
77 } else {
78 panel_str = idHelperSvc.toString(panelId) + " - INVALID ID";
79 }
80}
81
82//========================================================================================================
86
87//========================================================================================================
88std::string RpcPanel::getElementStr() const {
89 std::ostringstream ele_key;
90
91 ele_key << stationName << "_" << stationEta << "_" << stationPhi << "_"
92 << doubletR << "_" << doubletZ;
93 return ele_key.str();
94}
95
96//========================================================================================================
97bool RpcPanel::operator==(const RpcPanel &rhs) const {
98 if (this->stationName == rhs.stationName &&
99 this->stationPhi == rhs.stationPhi &&
100 this->stationEta == rhs.stationEta && this->doubletR == rhs.doubletR &&
101 this->doubletPhi == rhs.doubletPhi && this->doubletZ == rhs.doubletZ &&
102 this->gasGap == rhs.gasGap && this->measPhi == rhs.measPhi) {
103 return true;
104 } else {
105 return false;
106 }
107}
108
109//========================================================================================================
110bool RpcPanel::operator<(const RpcPanel &rhs) const {
111 return this->panel_index < rhs.panel_index;
112}
113
114//========================================================================================================
115bool RpcPanel::operator>(const RpcPanel &rhs) const {
116 return this->panel_index > rhs.panel_index;
117}
118
119//========================================================================================================
120std::string RpcPanel::getOnlineConvention() const {
121 const std::string mystring_name = stationNameStr;
122 const std::string mystring_eta = std::to_string(std::abs(stationEta));
123 const std::string mystring_dPhi = "DP" + std::to_string(doubletPhi) + ".";
124 const std::string mystring_gap = "Ly" + std::to_string(gasGap - 1) + ".";
125 const std::string mystring_dZ = "DZ" + std::to_string(doubletZ);
126
127 std::string mystring_phi;
128 std::string mystring_side;
129 std::string mystring_PICO; // pivot or confirm
130 std::string mystring_ETAPHI;
131
132 int myphi_part = (2 * stationPhi) - 1;
133
134 if (stationEta > 0) {
135 mystring_side = "A";
136 } else if (stationEta < 0) {
137 mystring_side = "C";
138 } else {
139 mystring_side = "N";
140 } // N: NULL, "ETA NOT DEFINED"
141
142 if (stationName == 3 || stationName == 5 || stationName == 8 ||
143 stationName == 9 || stationName == 10) {
144 myphi_part += 1;
145 }
146
147 if (myphi_part < 10) {
148 mystring_phi = "0" + std::to_string(myphi_part);
149 } else {
150 mystring_phi = std::to_string(myphi_part);
151 }
152
153 if ((stationNameStr == "BMS" && doubletR == 1) ||
154 (stationNameStr == "BML" && doubletR == 1) ||
155 (stationNameStr == "BMF" && doubletR == 1)) {
156 mystring_PICO = ".CO.";
157 } else if ((stationNameStr == "BMS" && doubletR == 2) ||
158 (stationNameStr == "BML" && doubletR == 2) ||
159 (stationNameStr == "BMF" && doubletR == 2)) {
160 mystring_PICO = ".PI.";
161 } else if ((stationNameStr == "BOS" && doubletR == 1) ||
162 (stationNameStr == "BOL" && doubletR == 1) ||
163 (stationNameStr == "BOF" && doubletR == 1) ||
164 (stationNameStr == "BOG" && doubletR == 1)) {
165 mystring_PICO = ".CO.";
166 } else {
167 // mystring_PICO = ".PIVOT_CONFIRMED_NOT_DEFINED.";
168 mystring_PICO =
169 ".NU" + std::to_string(doubletR) + "."; // NU: NULL, Not defined
170 }
171
172 if (measPhi == 0) {
173 mystring_ETAPHI = "ETA.";
174 } else {
175 mystring_ETAPHI = "PHI.";
176 }
177
178 return mystring_name + mystring_eta + mystring_side + mystring_phi +
179 mystring_PICO + mystring_dPhi + mystring_gap + mystring_dZ +
180 mystring_ETAPHI;
181}
182
183//========================================================================================================
184std::pair<int, int> RpcPanel::getSectorLayer() const {
185 // {'2':'BML', '3':'BMS', '4':'BOL', '5':'BOS', '8':'BMF' , '9':'BOF',
186 // '10':'BOG', '53':'BME'}
187 int sector = 0;
188 int layer = 0;
189
190 if (stationName == 2 || stationName == 53) {
191 sector = stationPhi * 2 - 1;
192 layer = (doubletR - 1) * 2 + gasGap;
193 } else if (stationName == 4) {
194 sector = stationPhi * 2 - 1;
195 layer = 4 + (doubletR - 1) * 2 + gasGap;
196 } else if (stationName == 3 || stationName == 8) {
197 sector = stationPhi * 2;
198 layer = (doubletR - 1) * 2 + gasGap;
199 } else if (stationName == 5 || stationName == 9 || stationName == 10) {
200 sector = stationPhi * 2;
201 layer = 4 + (doubletR - 1) * 2 + gasGap;
202 }
203 if (stationEta < 0)
204 sector = -1 * sector;
205
206 std::pair<int, int> sec_layer(sector, layer);
207
208 return sec_layer;
209}
210
211//========================================================================================================
212// strcuct ExResult
213//========================================================================================================
215 const Trk::PropDirection direction_)
216 : gasgap_id(gasgap_id_), direction(direction_) {}
217
218//========================================================================================================
219// GasGapData
220//========================================================================================================
222 const MuonGM::RpcReadoutElement *readoutEl_,
223 const int doubletZ_, const int doubletPhi_,
224 const unsigned gasgap_)
225 : readoutEl(readoutEl_),
226 doubletR(readoutEl_->getDoubletR()),
227 doubletZ(doubletZ_),
228 doubletPhi(doubletPhi_),
229 gasgap(gasgap_) {
230 const RpcIdHelper &rpcIdHelper = idHelperSvc.rpcIdHelper();
231 const Identifier readEl_id = readoutEl->identify();
232
233 stationName = rpcIdHelper.stationName(readEl_id);
234 stationEta = rpcIdHelper.stationEta(readEl_id);
235 stationPhi = rpcIdHelper.stationPhi(readEl_id);
236
237 //
238 // Get Identifier for this gas gap
239 //
242
243 if (gapid_valid) {
244 gapid_str = idHelperSvc.toString(gapid);
245 } else {
246 gapid_str = idHelperSvc.toString(gapid) + " - INVALID ID";
247 }
248
249 NetaStrips = readoutEl->NetaStrips();
250 NphiStrips = readoutEl->NphiStrips();
251
252 // Determine gas gap span in eta
253 const Identifier etaIdFirst{rpcIdHelper.channelID(gapid, doubletZ, doubletPhi, gasgap, false, 1)};
254 const Amg::Vector3D glbPos_etaStrip_1st = readoutEl->stripPos(etaIdFirst);
255 double eta_etaStrip_1st = glbPos_etaStrip_1st.eta();
256
257 const Identifier etaIdLast{rpcIdHelper.channelID(gapid, doubletZ, doubletPhi, gasgap, false, NetaStrips)};
258 const Amg::Vector3D glbPos_etaStrip_last = readoutEl->stripPos(etaIdLast);
259 double eta_etaStrip_last = glbPos_etaStrip_last.eta();
260
261 if (eta_etaStrip_1st < eta_etaStrip_last) {
262 minStripEta = eta_etaStrip_1st;
263 maxStripEta = eta_etaStrip_last;
264 } else {
265 minStripEta = eta_etaStrip_last;
266 maxStripEta = eta_etaStrip_1st;
267 }
268
269 //
270 // Determine gas gap span in phi
271 //
272 const Identifier phiIdFirst{rpcIdHelper.channelID(gapid, doubletZ, doubletPhi, gasgap, true, 1)};
273 const Amg::Vector3D glbPos_phiStrip_1st = readoutEl->stripPos(phiIdFirst);
274 double phi_phiStrip_1st = glbPos_phiStrip_1st.phi();
275
276 const Identifier phiIdLast{rpcIdHelper.channelID(gapid, doubletZ, doubletPhi, gasgap, true, NphiStrips)};
277 const Amg::Vector3D glbPos_phiStrip_last = readoutEl->stripPos(phiIdLast);
278 double phi_phiStrip_last = glbPos_phiStrip_last.phi();
279
280 if (phi_phiStrip_1st < phi_phiStrip_last) {
281 minStripPhi = phi_phiStrip_1st;
282 maxStripPhi = phi_phiStrip_last;
283 } else {
284 minStripPhi = phi_phiStrip_last;
285 maxStripPhi = phi_phiStrip_1st;
286 }
287}
288
289//========================================================================================================
291 ExResult &result, const xAOD::TrackParticle &track) const {
292 /*
293 This function:
294 - computes minum DR distance between track and RpcReadoutElement
295 - do this before expensive call to extrapolator
296 */
297
298 const double deta1 = std::abs(minStripEta - track.eta());
299 const double deta2 = std::abs(maxStripEta - track.eta());
300
301 const double dphi1 =
302 std::abs(xAOD::P4Helpers::deltaPhi(minStripPhi, track.phi()));
303 const double dphi2 =
304 std::abs(xAOD::P4Helpers::deltaPhi(maxStripPhi, track.phi()));
305
306 result.minTrackGasGapDEta = std::min<double>(deta1, deta2);
307 result.minTrackGasGapDPhi = std::min<double>(dphi1, dphi2);
308
309 //
310 // Now check if track position is between min and max - if true, set
311 // distance to zero
312 //
313 if (minStripEta <= track.eta() && track.eta() <= maxStripEta) {
314 result.minTrackGasGapDEta = 0.0;
315 }
316 if (minStripPhi <= track.phi() && track.phi() <= maxStripPhi) {
317 result.minTrackGasGapDPhi = 0.0;
318 }
319
320 //
321 // Compute min DR to this gas gap
322 //
323 result.minTrackGasGapDR =
324 std::sqrt(result.minTrackGasGapDEta * result.minTrackGasGapDEta +
325 result.minTrackGasGapDPhi * result.minTrackGasGapDPhi);
326}
327
328//========================================================================================================
330 ExResult &result, const Trk::TrackParameters *trackParam) const {
331 /*
332 This function:
333 - computes minum DR distance between track and RpcReadoutElement
334 - do this before expensive call to extrapolator
335 */
336
337 const double deta1 = std::abs(minStripEta - trackParam->position().eta());
338 const double deta2 = std::abs(maxStripEta - trackParam->position().eta());
339
340 const double dphi1 = std::abs(
341 xAOD::P4Helpers::deltaPhi(minStripPhi, trackParam->position().phi()));
342 const double dphi2 = std::abs(
343 xAOD::P4Helpers::deltaPhi(maxStripPhi, trackParam->position().phi()));
344
345 result.minTrackGasGapDEta = std::min<double>(deta1, deta2);
346 result.minTrackGasGapDPhi = std::min<double>(dphi1, dphi2);
347
348 //
349 // Now check if trackParam->position() position is between min and max - if
350 // true, set distance to zero
351 //
352 if (minStripEta <= trackParam->position().eta() &&
353 trackParam->position().eta() <= maxStripEta) {
354 result.minTrackGasGapDEta = 0.0;
355 }
356 if (minStripPhi <= trackParam->position().phi() &&
357 trackParam->position().phi() <= maxStripPhi) {
358 result.minTrackGasGapDPhi = 0.0;
359 }
360
361 //
362 // Compute min DR to this gas gap
363 //
364 result.minTrackGasGapDR =
365 std::sqrt(result.minTrackGasGapDEta * result.minTrackGasGapDEta +
366 result.minTrackGasGapDPhi * result.minTrackGasGapDPhi);
367}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
An RpcReadoutElement corresponds to a single RPC module; therefore typicaly a barrel muon station con...
int stationEta(const Identifier &id) const
int stationPhi(const Identifier &id) const
int stationName(const Identifier &id) const
const std::string & stationNameString(const int &index) const
Interface for Helper service that creates muon Identifiers and can be used to print Identifiers.
virtual const RpcIdHelper & rpcIdHelper() const =0
access to RpcIdHelper
virtual std::string toString(const Identifier &id) const =0
print all fields to string
Identifier channelID(int stationName, int stationEta, int stationPhi, int doubletR, int doubletZ, int doubletPhi, int gasGap, int measuresPhi, int strip) const
int gasGap(const Identifier &id) const override
get the hashes
Identifier panelID(const Identifier &padID, int gasGap, int measuresPhi) const
Identifier gapID(const Identifier &padID, int gasGap) const
int doubletPhi(const Identifier &id) const
int doubletR(const Identifier &id) const
bool measuresPhi(const Identifier &id) const override
int doubletZ(const Identifier &id) const
const Amg::Vector3D & position() const
Access method for the position.
Eigen::Matrix< double, 3, 1 > Vector3D
PropDirection
PropDirection, enum for direction of the propagation.
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
TrackParticle_v1 TrackParticle
Reference the current persistent version:
ExResult(const Identifier gasgap_id_, const Trk::PropDirection direction_)
const Trk::PropDirection direction
Definition RPCDQUtils.h:86
const Identifier gasgap_id
Definition RPCDQUtils.h:85
double minStripEta
Definition RPCDQUtils.h:124
Identifier gapid
Definition RPCDQUtils.h:110
void computeTrackDistanceToGasGap(ExResult &result, const xAOD::TrackParticle &track) const
unsigned NetaStrips
Definition RPCDQUtils.h:122
bool gapid_valid
Definition RPCDQUtils.h:111
const MuonGM::RpcReadoutElement * readoutEl
Definition RPCDQUtils.h:108
int stationName
Definition RPCDQUtils.h:114
double maxStripEta
Definition RPCDQUtils.h:125
unsigned NphiStrips
Definition RPCDQUtils.h:123
std::string gapid_str
Definition RPCDQUtils.h:112
double maxStripPhi
Definition RPCDQUtils.h:127
GasGapData(const Muon::IMuonIdHelperSvc &idHelperSvc, const MuonGM::RpcReadoutElement *readoutEl_, const int doubletZ_, const int doubletPhi_, const unsigned gasgap_)
double minStripPhi
Definition RPCDQUtils.h:126
unsigned gasgap
Definition RPCDQUtils.h:120
RpcPanel()=default
bool panel_valid
Definition RPCDQUtils.h:30
void SetPanelIndex(int index)
int stationEta
Definition RPCDQUtils.h:39
int doubletZ
Definition RPCDQUtils.h:42
bool operator==(const RpcPanel &rhs) const
std::string panel_name
Definition RPCDQUtils.h:32
const MuonGM::RpcReadoutElement * readoutEl
Definition RPCDQUtils.h:27
Identifier panelId
Definition RPCDQUtils.h:29
int gasGap
Definition RPCDQUtils.h:44
bool operator>(const RpcPanel &rhs) const
int doubletR
Definition RPCDQUtils.h:41
int panel_index
Definition RPCDQUtils.h:46
int measPhi
Definition RPCDQUtils.h:45
std::string stationNameStr
Definition RPCDQUtils.h:33
int stationPhi
Definition RPCDQUtils.h:40
bool operator<(const RpcPanel &rhs) const
int doubletPhi
Definition RPCDQUtils.h:43
void FillRpcId(Identifier id, const RpcIdHelper &rpcIdHelper)
int stationName
Definition RPCDQUtils.h:38
std::string panel_str
Definition RPCDQUtils.h:31
std::string getOnlineConvention() const
std::string getElementStr() const
std::pair< int, int > getSectorLayer() const