ATLAS Offline Software
Loading...
Searching...
No Matches
RpcRoadDefiner.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#include <cmath>
6#include <functional>
7#include <type_traits>
8
9#include "RpcRoadDefiner.h"
12
13// --------------------------------------------------------------------------------
14// --------------------------------------------------------------------------------
15
17{
18 ATH_CHECK(m_idHelperSvc.retrieve());
19
20 ATH_CHECK(m_regionSelector.retrieve());
21 ATH_MSG_DEBUG("Retrieved the RegionSelector tool ");
22
23 return StatusCode::SUCCESS;
24}
25
26// --------------------------------------------------------------------------------
27// --------------------------------------------------------------------------------
28
29StatusCode TrigL2MuonSA::RpcRoadDefiner::defineRoad(const EventContext& ctx,
30 const xAOD::MuonRoI* p_roi,
31 const bool insideOut,
32 TrigL2MuonSA::MuonRoad& muonRoad,
33 const TrigL2MuonSA::RpcLayerHits& rpcLayerHits,
34 const ToolHandle<RpcPatFinder>* rpcPatFinder,
35 TrigL2MuonSA::RpcFitResult& rpcFitResult,
36 const double roiEtaMinLow,
37 const double roiEtaMaxLow,
38 const double roiEtaMinHigh,
39 const double roiEtaMaxHigh) const
40{
41 const double ZERO_LIMIT = 1e-5;
42
43 if (m_use_rpc && !insideOut) {
44 std::array<std::reference_wrapper<double>, 3> aw {
45 rpcFitResult.slope_inner, rpcFitResult.slope_middle, rpcFitResult.slope_outer
46 };
47 std::array<std::reference_wrapper<double>, 3> bw {
48 rpcFitResult.offset_inner, rpcFitResult.offset_middle, rpcFitResult.offset_outer
49 };
50 rpcFitResult.isSuccess = (*rpcPatFinder)->findPatternEta(aw, bw, rpcLayerHits);
51
52 for(int i=0;i<3;i++){
53 if(std::abs(aw[i].get()) <= ZERO_LIMIT) rpcFitResult.isSuccess = false;
54 }
55
56 double phi_middle;
57 double phi_outer;
58 if ( (*rpcPatFinder)->findPatternPhi(phi_middle, phi_outer, rpcLayerHits)) {
59 rpcFitResult.phi = phi_middle;
60 rpcFitResult.phi_middle = phi_middle;
61 rpcFitResult.phi_outer = phi_outer;
62 } else {
63 rpcFitResult.phi = p_roi->phi();
64 }
65 ATH_MSG_DEBUG("RpcPatFinder: Found: " << rpcFitResult.isSuccess << " Slopes: " << aw[0] << "," << aw[1] << "," << aw[2] << " Offsets: " << bw[0] << "," << bw[1] << "," << bw[2]);
66 } else {
67 ATH_MSG_DEBUG("Skip rpcPatFinder");
68 }
69
70 muonRoad.isEndcap = false;
71 if(!insideOut) {
72 muonRoad.phiMiddle = rpcFitResult.phi;
73 } else {
74 muonRoad.phiMiddle = muonRoad.extFtfMiddlePhi;
75 rpcFitResult.phi = muonRoad.extFtfMiddlePhi;
76 rpcFitResult.phi_middle = muonRoad.extFtfMiddlePhi;
77 rpcFitResult.phi_outer = muonRoad.extFtfMiddlePhi;
78 }
79 muonRoad.phiRoI = p_roi->phi();
80 muonRoad.side = (p_roi->phi()<0.)? 0 : 1;
81 muonRoad.LargeSmall = ((p_roi->getSectorID() + 1)/2 )%2;
82
83 const int PhysicsSector = ((p_roi->getSectorID() + 1)/4 )%8 + 1;
84
85 int special = 0;
86 if (muonRoad.LargeSmall == 0 && (PhysicsSector == 6 || PhysicsSector == 8 )) special = 1; // BIM BIR
87 if (muonRoad.LargeSmall == 1 && (PhysicsSector == 6 || PhysicsSector == 7 )) special = 1; //feets
88 muonRoad.Special = special;
89
90 auto fillAllLayersWith = [&muonRoad](const int& station, const double& value) -> void {
91 std::fill(std::begin(muonRoad.rWidth[station]), std::end(muonRoad.rWidth[station]), value);
92 };
93
94 if (!rpcFitResult.isSuccess) {
95 if (!insideOut) {
96 fillAllLayersWith(0, 500); // Inner
97 fillAllLayersWith(1, 650); // Middle
98 fillAllLayersWith(2, 800); // Outer
99 fillAllLayersWith(3, 500); // EndcapInner
100 fillAllLayersWith(9, 650); // BME
101 fillAllLayersWith(10, 650); // BMG
102 } else {
103 fillAllLayersWith(0, 250); // Inner
104 fillAllLayersWith(1, 400); // Middle
105 fillAllLayersWith(2, 600); // Outer
106 fillAllLayersWith(3, 300); // EndcapInner
107 fillAllLayersWith(9, 400); // BME
108 fillAllLayersWith(10, 400); // BMG
109 }
110 } else {
111 fillAllLayersWith(0, 400); // Inner
112 fillAllLayersWith(1, 200); // Middle
113 fillAllLayersWith(2, 400); // Outer
114 fillAllLayersWith(3, 400); // EndcapInner
115 fillAllLayersWith(9, m_rWidth_RPC_Failed); // BME
116 fillAllLayersWith(10, m_rWidth_RPC_Failed); // BMG
117 }
118
119 std::vector<IdentifierHash> mdtHashList;
120
121 // get sector_trigger and sector_overlap by using the region selector
122 IdContext context = m_idHelperSvc->mdtIdHelper().module_context();
123
124 double etaMin = p_roi->eta()-.02;
125 double etaMax = p_roi->eta()+.02;
126 double phiMin = muonRoad.phiMiddle-.01;
127 double phiMax = muonRoad.phiMiddle+.01;
128 if(phiMax > M_PI) phiMax -= M_PI*2.;
129 if(phiMin < -M_PI) phiMin += M_PI*2.;
130
131 auto roi = std::make_unique<TrigRoiDescriptor>( p_roi->eta(), etaMin, etaMax, muonRoad.phiMiddle, phiMin, phiMax );
132
133 if (roi) m_regionSelector->lookup(ctx)->HashIDList( *roi, mdtHashList);
134 else {
135 TrigRoiDescriptor fullscan_roi( true );
136 m_regionSelector->lookup(ctx)->HashIDList(fullscan_roi, mdtHashList);
137 }
138
139 int &sector_trigger {muonRoad.MDT_sector_trigger}, &sector_overlap {muonRoad.MDT_sector_overlap};
140 sector_trigger = (PhysicsSector - 1)*2 + muonRoad.LargeSmall;
141 sector_overlap = 99;
142
143 for( const IdentifierHash& hash : mdtHashList){
144
145 Identifier id;
146 if( m_idHelperSvc->mdtIdHelper().get_id(hash, id, &context) !=0 ) ATH_MSG_ERROR("problem converting hash list to id");
147
148 muonRoad.stationList.push_back(id);
149 const int stationPhi = m_idHelperSvc->mdtIdHelper().stationPhi(id);
150 const std::string name = m_idHelperSvc->mdtIdHelper().stationNameString(m_idHelperSvc->mdtIdHelper().stationName(id));
151
152 if ( name[1]=='M' && name[2]=='E' ) continue;//exclude BME
153 if ( name[1]=='M' && name[2]=='G' ) continue;//exclude BMG
154
155 const int LargeSmall = (name[2]=='S' || name[2]=='F' || name[2]=='G' ) ? 1 : 0;
156 const int sector = (stationPhi-1)*2 + LargeSmall;
157
158 if (sector != sector_trigger) sector_overlap = sector;
159 if (sector != sector_trigger and sector_overlap != 99 and sector != sector_overlap) ATH_MSG_ERROR("Multiple sector overlaps not expected");
160 }
161
163 auto fillAllSectorsWith = [&muonRoad](const int& station, const double& aw, const double& bw) -> void {
164 std::fill(std::begin(muonRoad.aw[station]), std::end(muonRoad.aw[station]), aw);
165 std::fill(std::begin(muonRoad.bw[station]), std::end(muonRoad.bw[station]), bw);
166 };
167
168 if (!insideOut) {
169 if (rpcFitResult.isSuccess) { // MOD this code is very suspicious
170 fillAllSectorsWith(0, rpcFitResult.slope_inner, rpcFitResult.offset_inner); //Barrel inner
171 fillAllSectorsWith(1, rpcFitResult.slope_middle, rpcFitResult.offset_middle); //Barrel middle
172 fillAllSectorsWith(2, rpcFitResult.slope_outer, rpcFitResult.offset_outer); //Barrel outer
173 fillAllSectorsWith(3, rpcFitResult.slope_inner, rpcFitResult.offset_inner); //Endcap inner
174 fillAllSectorsWith(9, rpcFitResult.slope_middle, rpcFitResult.offset_middle); //BME
175 fillAllSectorsWith(10, rpcFitResult.slope_middle, rpcFitResult.offset_middle); //BMG
176
177 } else {
178 auto compute_aw = [&ZERO_LIMIT](double etaMin, double etaMax) -> double {
179 const double eta = 0.5 * (etaMin + etaMax);
180 if (std::abs(eta) < ZERO_LIMIT) return 0.;
181 const double theta = 2.* std::atan(std::exp(-std::abs(eta)));
182 return std::tan(theta) * (eta / std::abs(eta)); // preserves sign
183 };
184
185 const double awLow = compute_aw(roiEtaMinLow, roiEtaMaxLow);
186 const double awHigh = compute_aw(roiEtaMinHigh, roiEtaMaxHigh);
187
188 fillAllSectorsWith(0, awLow, 0.); //Barrel inner
189 fillAllSectorsWith(1, awLow, 0.); //Barrel middle
190 fillAllSectorsWith(2, awHigh, 0.); //Barrel outer
191 fillAllSectorsWith(3, awLow, 0.); //Endcap inner
192 fillAllSectorsWith(9, awLow, 0.); //BME
193 fillAllSectorsWith(10, awLow, 0.); //BMG
194 }
195 } else {
196
197 ATH_MSG_DEBUG("Use aw_ftf and bw_ftf as aw and bw");
198 fillAllSectorsWith(0, muonRoad.aw_ftf[0][0], muonRoad.bw_ftf[0][0]); //Barrel inner
199 fillAllSectorsWith(1, muonRoad.aw_ftf[1][0], muonRoad.bw_ftf[1][0]); //Barrel middle
200 fillAllSectorsWith(2, muonRoad.aw_ftf[2][0], muonRoad.bw_ftf[2][0]); //Barrel outer
201 fillAllSectorsWith(3, muonRoad.aw_ftf[3][0], muonRoad.bw_ftf[3][0]); //Endcap inner
202 fillAllSectorsWith(9, muonRoad.aw_ftf[9][0], muonRoad.bw_ftf[9][0]); //BME
203 fillAllSectorsWith(10, muonRoad.aw_ftf[10][0], muonRoad.bw_ftf[10][0]); //BMG
204 }
205
206 ATH_MSG_DEBUG("muonRoad.phiMiddle: " << muonRoad.phiMiddle);
207
208 return StatusCode::SUCCESS;
209}
210
211// --------------------------------------------------------------------------------
212// --------------------------------------------------------------------------------
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
const float ZERO_LIMIT
This class saves the "context" of an expanded identifier (ExpandedIdentifier) for compact or hash ver...
Definition IdContext.h:26
This is a "hash" representation of an Identifier.
std::vector< Identifier > stationList
Definition MuonRoad.h:102
double aw[N_STATION][N_SECTOR]
Definition MuonRoad.h:83
double bw_ftf[N_STATION][N_SECTOR]
Definition MuonRoad.h:95
double aw_ftf[N_STATION][N_SECTOR]
Definition MuonRoad.h:94
double bw[N_STATION][N_SECTOR]
Definition MuonRoad.h:84
double rWidth[N_STATION][N_LAYER]
Definition MuonRoad.h:86
StatusCode defineRoad(const EventContext &ctx, const xAOD::MuonRoI *p_roi, const bool insideOut, TrigL2MuonSA::MuonRoad &muonRoad, const TrigL2MuonSA::RpcLayerHits &rpcLayerHits, const ToolHandle< RpcPatFinder > *rpcPatFinder, TrigL2MuonSA::RpcFitResult &rpcFitResult, const double roiEtaMinLow, const double roiEtaMaxLow, const double roiEtaMinHigh, const double roiEtaMaxHigh) const
ToolHandle< IRegSelTool > m_regionSelector
virtual StatusCode initialize() override
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
nope - should be used for standalone also, perhaps need to protect the class def bits ifndef XAOD_ANA...
float eta() const
The pseudorapidity ( ) of the muon candidate.
float phi() const
The azimuthal angle ( ) of the muon candidate.
int getSectorID() const
Get the sector ID number.
T * get(TKey *tobj)
get a TObject* from a TKey* (why can't a TObject be a TKey?)
Definition hcg.cxx:130
MuonRoI_v1 MuonRoI
Definition MuonRoI.h:15