ATLAS Offline Software
Loading...
Searching...
No Matches
TrigL2MuonSA::RpcRoadDefiner Class Reference

#include <RpcRoadDefiner.h>

Inheritance diagram for TrigL2MuonSA::RpcRoadDefiner:
Collaboration diagram for TrigL2MuonSA::RpcRoadDefiner:

Public Member Functions

virtual StatusCode initialize () override
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
void setRoadWidthForFailure (double rWidth_RPC_Failed)
void setRpcGeometry (bool use_rpc)
 AthAlgTool (const std::string &type, const std::string &name, const IInterface *parent)
 Constructor with parameters:
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

float f (float x, float c0, float c1, float c2, float c3) const
float fp (float x, float c33, float c22, float c1) const
void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

double m_rWidth_RPC_Failed {0}
bool m_use_rpc {true}
ToolHandle< IRegSelToolm_regionSelector {this, "RegionSelectionTool", "RegSelTool/RegSelTool_MDT", "MDT Region Selector Tool"}
ServiceHandle< Muon::IMuonIdHelperSvcm_idHelperSvc {this, "MuonIdHelperSvc", "Muon::MuonIdHelperSvc/MuonIdHelperSvc"}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 30 of file RpcRoadDefiner.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Member Function Documentation

◆ AthAlgTool()

AthAlgTool::AthAlgTool ( const std::string & type,
const std::string & name,
const IInterface * parent )

Constructor with parameters:

Definition at line 31 of file AthAlgTool.cxx.

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ defineRoad()

StatusCode TrigL2MuonSA::RpcRoadDefiner::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

Copy a,b coefficients to muonRoad

Definition at line 29 of file RpcRoadDefiner.cxx.

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}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar theta() const
theta method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Current TrigRoiDescriptor
const float ZERO_LIMIT
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
ToolHandle< IRegSelTool > m_regionSelector
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
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
constexpr uint8_t stationPhi
station Phi 1 to 8

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ f()

float TrigL2MuonSA::RpcRoadDefiner::f ( float x,
float c0,
float c1,
float c2,
float c3 ) const
inlineprotected

Definition at line 68 of file RpcRoadDefiner.h.

69 {
70 return c0 + x * (c1 + x * (c2 + x * c3)); // faster
71 }
#define x

◆ fp()

float TrigL2MuonSA::RpcRoadDefiner::fp ( float x,
float c33,
float c22,
float c1 ) const
inlineprotected

Definition at line 75 of file RpcRoadDefiner.h.

76 {
77 return c1 + x * (c22 + x * c33); // faster
78 }

◆ initialize()

StatusCode TrigL2MuonSA::RpcRoadDefiner::initialize ( )
overridevirtual

Definition at line 16 of file RpcRoadDefiner.cxx.

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}
#define ATH_CHECK
Evaluate an expression and check for errors.

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ setRoadWidthForFailure()

void TrigL2MuonSA::RpcRoadDefiner::setRoadWidthForFailure ( double rWidth_RPC_Failed)
inline

Definition at line 50 of file RpcRoadDefiner.h.

50{ m_rWidth_RPC_Failed = rWidth_RPC_Failed; };

◆ setRpcGeometry()

void TrigL2MuonSA::RpcRoadDefiner::setRpcGeometry ( bool use_rpc)
inline

Definition at line 51 of file RpcRoadDefiner.h.

51{ m_use_rpc = use_rpc; };

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_idHelperSvc

ServiceHandle<Muon::IMuonIdHelperSvc> TrigL2MuonSA::RpcRoadDefiner::m_idHelperSvc {this, "MuonIdHelperSvc", "Muon::MuonIdHelperSvc/MuonIdHelperSvc"}
private

Definition at line 62 of file RpcRoadDefiner.h.

62{this, "MuonIdHelperSvc", "Muon::MuonIdHelperSvc/MuonIdHelperSvc"};

◆ m_regionSelector

ToolHandle<IRegSelTool> TrigL2MuonSA::RpcRoadDefiner::m_regionSelector {this, "RegionSelectionTool", "RegSelTool/RegSelTool_MDT", "MDT Region Selector Tool"}
private

Definition at line 61 of file RpcRoadDefiner.h.

61{this, "RegionSelectionTool", "RegSelTool/RegSelTool_MDT", "MDT Region Selector Tool"};

◆ m_rWidth_RPC_Failed

double TrigL2MuonSA::RpcRoadDefiner::m_rWidth_RPC_Failed {0}
private

Definition at line 58 of file RpcRoadDefiner.h.

58{0};

◆ m_use_rpc

bool TrigL2MuonSA::RpcRoadDefiner::m_use_rpc {true}
private

Definition at line 59 of file RpcRoadDefiner.h.

59{true};

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


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