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

#include <ClusterRoadDefiner.h>

Inheritance diagram for TrigL2MuonSA::ClusterRoadDefiner:
Collaboration diagram for TrigL2MuonSA::ClusterRoadDefiner:

Public Member Functions

 ClusterRoadDefiner (const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode initialize () override
StatusCode defineRoad (const EventContext &ctx, const xAOD::MuonRoI *p_roi, std::vector< TrigL2MuonSA::MuonRoad > &clusterRoad, TrigL2MuonSA::RpcLayerClusters &rpcLayerClusters, const ToolHandle< ClusterPatFinder > *clusterPatFinder, std::vector< TrigL2MuonSA::RpcFitResult > &clusterFitResults, double roiEtaMinLow, double roiEtaMaxLow, double roiEtaMinHigh, double roiEtaMaxHigh) const
void setRoadWidthForFailure (double rWidth_RPC_Failed)
void setRpcGeometry (bool use_rpc)
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

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::MuonIdHleperSvc/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 33 of file ClusterRoadDefiner.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ ClusterRoadDefiner()

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

Definition at line 13 of file ClusterRoadDefiner.cxx.

15 :
16 AthAlgTool(type, name, parent)
17{
18}
AthAlgTool()
Default constructor:

Member Function Documentation

◆ 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::ClusterRoadDefiner::defineRoad ( const EventContext & ctx,
const xAOD::MuonRoI * p_roi,
std::vector< TrigL2MuonSA::MuonRoad > & clusterRoad,
TrigL2MuonSA::RpcLayerClusters & rpcLayerClusters,
const ToolHandle< ClusterPatFinder > * clusterPatFinder,
std::vector< TrigL2MuonSA::RpcFitResult > & clusterFitResults,
double roiEtaMinLow,
double roiEtaMaxLow,
double roiEtaMinHigh,
double roiEtaMaxHigh ) const

Definition at line 35 of file ClusterRoadDefiner.cxx.

45{
46 const double ZERO_LIMIT = 1e-5;
47 const int N_STATION = 3; //0:inner, 1:middle, 2:outer
48
49 const int N_LAYER = 5; // 0: inner, 1: middle, 2: outer 4: BME 5: BMG
50 const int N_SECTOR = 2; // 0: normal, 1:overlap
51
52 if (m_use_rpc) {
53 std::vector< std::vector< double > > aw, bw;
54 std::vector< double > clearRoad;
55 clearRoad.clear();
56 aw.assign(N_STATION, clearRoad);
57 bw.assign(N_STATION, clearRoad);
58 ATH_MSG_DEBUG("start searching eta pattern");
59 if ( (*clusterPatFinder)->findPatternEta(aw, bw, rpcLayerClusters) ) {
60 for(unsigned int iClus = 0; iClus < aw[1].size(); iClus++){
61 TrigL2MuonSA::RpcFitResult clusFitResult;
62 clusFitResult.Clear();
63 clusFitResult.isSuccess = true;
64 clusFitResult.offset_inner = (bw[0]).at(iClus);
65 clusFitResult.offset_middle = (bw[1]).at(iClus);
66 clusFitResult.offset_outer = (bw[2]).at(iClus);
67 clusFitResult.slope_inner = 1.0/(aw[0]).at(iClus);
68 clusFitResult.slope_middle = 1.0/(aw[1]).at(iClus);
69 clusFitResult.slope_outer = 1.0/(aw[2]).at(iClus);
70 ATH_MSG_DEBUG("==========================================================");
71 ATH_MSG_DEBUG("inner clusterRoad slope/offset = " << clusFitResult.slope_inner << "/" << clusFitResult.offset_inner);
72 ATH_MSG_DEBUG("middle clusterRoad slope/offset = " << clusFitResult.slope_middle << "/" << clusFitResult.offset_middle);
73 ATH_MSG_DEBUG("outer clusterRoad slope/offset = " << clusFitResult.slope_outer << "/" << clusFitResult.offset_outer);
74 for(int i=0;i<N_STATION;i++){
75 if(std::abs(1.0/aw[i].at(iClus)) <= ZERO_LIMIT) clusFitResult.isSuccess = false;
76 }
77 clusterFitResults.push_back(clusFitResult);
78 }
79 ATH_MSG_DEBUG("==========================================================");
80 } else {
81 TrigL2MuonSA::RpcFitResult clusFitResult;
82 clusFitResult.Clear();
83 clusFitResult.isSuccess = false;
84 clusterFitResults.push_back(clusFitResult);
85 }
86
87 ATH_MSG_DEBUG("start searching phi pattern");
88 std::vector<double> phi_middle, phi_outer;
89 phi_middle.clear(); phi_outer.clear();
90 if ( (*clusterPatFinder)->findPatternPhi(phi_middle, phi_outer, rpcLayerClusters) ){
91 double phi_middle_tot = 0;
92 double phi_outer_tot = 0;
93 size_t npatternPhi = phi_middle.size();
94 if(npatternPhi > 0){
95 for(double& phi : phi_middle){
96 phi_middle_tot = phi_middle_tot + phi;
97 }
98 for(double& phi : phi_outer){
99 phi_outer_tot = phi_outer_tot + phi;
100 }
101 double phi_middle_center = phi_middle_tot/npatternPhi;
102 double phi_outer_center = phi_outer_tot/npatternPhi;
103 ATH_MSG_DEBUG("center of phi middle/outer = " << phi_middle_center << "/" << phi_outer_center);
104 for(unsigned int iClus_fit = 0; iClus_fit < clusterFitResults.size(); iClus_fit++){
105 clusterFitResults.at(iClus_fit).phi = phi_middle_center;
106 clusterFitResults.at(iClus_fit).phi_middle = phi_middle_center;
107 clusterFitResults.at(iClus_fit).phi_outer = phi_outer_center;
108 }
109 } else {
110 for(unsigned int iClus_fit = 0; iClus_fit < clusterFitResults.size(); iClus_fit++){
111 clusterFitResults.at(iClus_fit).phi = p_roi->phi();
112 clusterFitResults.at(iClus_fit).phi_middle = p_roi->phi();
113 clusterFitResults.at(iClus_fit).phi_outer = p_roi->phi();
114 }
115 }
116 } else {
117 for(unsigned int iClus_fit = 0; iClus_fit < clusterFitResults.size(); iClus_fit++){
118 clusterFitResults.at(iClus_fit).phi = p_roi->phi();
119 }
120 }
121 }
122 if(clusterFitResults.empty()){
123 TrigL2MuonSA::RpcFitResult clusFitResult;
124 clusFitResult.Clear();
125 clusFitResult.isSuccess = false;
126 clusFitResult.phi = p_roi->phi();
127 clusterFitResults.push_back(clusFitResult);
128 }
129
130 ATH_MSG_DEBUG("stored cluster eta/phi pattern to clusterRoad");
131 for(unsigned int iClus_fit = 0; iClus_fit < clusterFitResults.size(); iClus_fit++){
132 TrigL2MuonSA::MuonRoad muonRoad;
133 muonRoad.Clear();
134 // RPC data is not available -> use RoI
135
136 muonRoad.isEndcap = false;
137 muonRoad.phiMiddle = clusterFitResults.at(iClus_fit).phi;
138 muonRoad.phiRoI = p_roi->phi();
139 muonRoad.side = (p_roi->phi()<0.)? 0 : 1;
140 muonRoad.LargeSmall = ((p_roi->getSectorID() + 1)/2 )%2;
141
142 int PhysicsSector = ((p_roi->getSectorID() + 1)/4 )%8 + 1;
143
144 int special = 0;
145 if (muonRoad.LargeSmall == 0 && (PhysicsSector == 6 || PhysicsSector == 8 ))
146 special = 1;
147 else if (muonRoad.LargeSmall == 1 && (PhysicsSector == 6 || PhysicsSector == 7 ))
148 special = 1;
149 muonRoad.Special = special;
150
151 for (int i_station=0; i_station<6; i_station++) {
152 for (int i_layer=0; i_layer<8; i_layer++) {
153 if(!clusterFitResults.at(iClus_fit).isSuccess) {
154 if (i_station==0) muonRoad.rWidth[i_station][i_layer] = 500;//for inner
155 else if (i_station==1) muonRoad.rWidth[i_station][i_layer] = 650;//for middle
156 else if (i_station==2) muonRoad.rWidth[i_station][i_layer] = 800;//for outer
157 else if (i_station==3) muonRoad.rWidth[i_station][i_layer] = 500;//EndcapInner
158 else if (i_station==4) muonRoad.rWidth[9][i_layer] = 650;//BME
159 else if (i_station==5) muonRoad.rWidth[10][i_layer] = 650;//BMG
160 else muonRoad.rWidth[i_station][i_layer] = m_rWidth_RPC_Failed;
161 }
162 else {
163 if (i_station==0) muonRoad.rWidth[i_station][i_layer] = 400;//for inner
164 else if (i_station==1) muonRoad.rWidth[i_station][i_layer] = 200;//for middle
165 else if (i_station==2) muonRoad.rWidth[i_station][i_layer] = 400;//for outer
166 else if (i_station==3) muonRoad.rWidth[i_station][i_layer] = 400;//EndcapInner
167 else if (i_station==4) muonRoad.rWidth[9][i_layer] = m_rWidth_RPC_Failed;//BME
168 else if (i_station==5) muonRoad.rWidth[10][i_layer] = m_rWidth_RPC_Failed;//BMG
169 else muonRoad.rWidth[i_station][i_layer] = m_rWidth_RPC_Failed;
170 }
171 }
172 }
173 int sector_trigger = 99;
174 int sector_overlap = 99;
175 std::vector<Identifier> stationList;
176 std::vector<IdentifierHash> mdtHashList;
177
178 // get sector_trigger and sector_overlap by using the region selector
179 IdContext context = m_idHelperSvc->mdtIdHelper().module_context();
180
181 double etaMin = p_roi->eta()-.02;
182 double etaMax = p_roi->eta()+.02;
183 double phiMin = muonRoad.phiMiddle-.01;
184 double phiMax = muonRoad.phiMiddle+.01;
185 if(phiMax > M_PI) phiMax -= M_PI*2.;
186 if(phiMin < M_PI*-1) phiMin += M_PI*2.;
187
188 TrigRoiDescriptor* roi = new TrigRoiDescriptor( p_roi->eta(), etaMin, etaMax, p_roi->phi(), phiMin, phiMax );
189
190 const IRoiDescriptor* iroi = static_cast<IRoiDescriptor*> (roi);
191
192 if (iroi) m_regionSelector->lookup(ctx)->HashIDList(*iroi, mdtHashList);
193 else {
194 TrigRoiDescriptor fullscan_roi( true );
195 m_regionSelector->lookup(ctx)->HashIDList(fullscan_roi, mdtHashList);
196 }
197
198 if(roi) delete roi;
199
200 for( const IdentifierHash& hash : mdtHashList){
201
202 Identifier id;
203 const int convert = m_idHelperSvc->mdtIdHelper().get_id(hash, id, &context);
204
205 if(convert!=0) ATH_MSG_ERROR("problem converting hash list to id");
206
207 muonRoad.stationList.push_back(id);
208 const int stationPhi = m_idHelperSvc->mdtIdHelper().stationPhi(id);
209 std::string name = m_idHelperSvc->mdtIdHelper().stationNameString(m_idHelperSvc->mdtIdHelper().stationName(id));
210
211 if ( name[1]=='M' && name[2]=='E' ) continue;//exclude BME
212 else if ( name[1]=='M' && name[2]=='G' ) continue;//exclude BMG
213
214 int LargeSmall = 0;
215 if(name[2]=='S' || name[2]=='F' || name[2]=='G' ) LargeSmall = 1;
216 int sector = (stationPhi-1)*2 + LargeSmall;
217 if(sector_trigger == 99)
218 sector_trigger = sector;
219 else if(sector_trigger != sector)
220 sector_overlap = sector;
221 }
222
223 int MDT_tr = (PhysicsSector - 1)*2 + muonRoad.LargeSmall;
224 if (MDT_tr == sector_overlap) {
225 sector_overlap = sector_trigger;
226 sector_trigger = MDT_tr;
227 }
228
229 muonRoad.MDT_sector_trigger = sector_trigger;
230 muonRoad.MDT_sector_overlap = sector_overlap;
231
232 if (clusterFitResults.at(iClus_fit).isSuccess) {
233 for (int i_sector=0; i_sector<N_SECTOR; i_sector++) {
234 muonRoad.aw[0][i_sector] = clusterFitResults.at(iClus_fit).slope_inner;
235 muonRoad.bw[0][i_sector] = clusterFitResults.at(iClus_fit).offset_inner;
236 muonRoad.aw[1][i_sector] = clusterFitResults.at(iClus_fit).slope_middle;
237 muonRoad.bw[1][i_sector] = clusterFitResults.at(iClus_fit).offset_middle;
238 muonRoad.aw[2][i_sector] = clusterFitResults.at(iClus_fit).slope_outer;
239 muonRoad.bw[2][i_sector] = clusterFitResults.at(iClus_fit).offset_outer;
240 muonRoad.aw[3][i_sector] = clusterFitResults.at(iClus_fit).slope_inner; // Endcap Inner
241 muonRoad.bw[3][i_sector] = clusterFitResults.at(iClus_fit).offset_inner;
242 muonRoad.aw[9][i_sector] = clusterFitResults.at(iClus_fit).slope_middle;//BME
243 muonRoad.bw[9][i_sector] = clusterFitResults.at(iClus_fit).offset_middle;
244 muonRoad.aw[10][i_sector] = clusterFitResults.at(iClus_fit).slope_middle;//BMG
245 muonRoad.bw[10][i_sector] = clusterFitResults.at(iClus_fit).offset_middle;
246 }
247 } else {
248 double roiEtaLow = (roiEtaMinLow + roiEtaMaxLow) * 0.5;
249 double roiEtaHigh = (roiEtaMinHigh + roiEtaMaxHigh) * 0.5;
250 double thetaLow = std::atan(std::exp(-std::abs(roiEtaLow)))*2.;
251 double thetaHigh = std::atan(std::exp(-std::abs(roiEtaHigh)))*2.;
252 double awLow = (std::abs(roiEtaLow) > ZERO_LIMIT)? std::tan(thetaLow)*(std::abs(roiEtaLow)/roiEtaLow): 0.;
253 double awHigh = (std::abs(roiEtaHigh) > ZERO_LIMIT)? std::tan(thetaHigh)*(std::abs(roiEtaHigh)/roiEtaHigh): 0.;
254
255 for (int i_station=0; i_station<N_LAYER; i_station++) {
256 for (int i_sector=0; i_sector<N_SECTOR; i_sector++) {
257 muonRoad.aw[i_station][i_sector] = awLow;
258 muonRoad.bw[i_station][i_sector] = 0;
259 if (i_station==2) muonRoad.aw[i_station][i_sector] = awHigh;
260 else if (i_station==3) muonRoad.aw[i_station][i_sector] = awLow; //EI
261 else if (i_station==4) muonRoad.aw[9][i_sector] = awLow; //BME
262 else if (i_station==5) muonRoad.aw[10][i_sector] = awLow; //BMG
263 }
264 }
265 }
266 clusterRoad.push_back(muonRoad);
267 }
268 ATH_MSG_DEBUG("finished ClusterRoadDefiner algorithm... leave");
269 return StatusCode::SUCCESS;
270}
#define M_PI
Scalar phi() const
phi method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Current TrigRoiDescriptor
const float ZERO_LIMIT
ToolHandle< IRegSelTool > m_regionSelector
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
std::vector< Identifier > stationList
Definition MuonRoad.h:102
double aw[N_STATION][N_SECTOR]
Definition MuonRoad.h:83
double bw[N_STATION][N_SECTOR]
Definition MuonRoad.h:84
double rWidth[N_STATION][N_LAYER]
Definition MuonRoad.h:86
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.
constexpr uint8_t stationPhi
station Phi 1 to 8
std::unique_ptr< MVAUtils::BDT > convert(TMVA::MethodBDT *bdt, bool isRegression=true, bool useYesNoLeaf=false)
constexpr int N_STATION
Definition MuonRoad.h:15
constexpr int N_LAYER
Definition MuonRoad.h:17
constexpr int N_SECTOR
Definition MuonRoad.h:16

◆ 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

◆ initialize()

StatusCode TrigL2MuonSA::ClusterRoadDefiner::initialize ( )
overridevirtual

Definition at line 23 of file ClusterRoadDefiner.cxx.

24{
25 ATH_CHECK(m_idHelperSvc.retrieve());
26 ATH_CHECK(m_regionSelector.retrieve());
27 ATH_MSG_DEBUG("Retrieved the RegionSelector tool ");
28
29 return StatusCode::SUCCESS;
30}
#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::ClusterRoadDefiner::setRoadWidthForFailure ( double rWidth_RPC_Failed)
inline

Definition at line 53 of file ClusterRoadDefiner.h.

53{ m_rWidth_RPC_Failed = rWidth_RPC_Failed; };

◆ setRpcGeometry()

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

Definition at line 54 of file ClusterRoadDefiner.h.

54{ 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::ClusterRoadDefiner::m_idHelperSvc {this, "MuonIdHelperSvc", "Muon::MuonIdHleperSvc/MuonIdHelperSvc"}
private

Definition at line 61 of file ClusterRoadDefiner.h.

61{this, "MuonIdHelperSvc", "Muon::MuonIdHleperSvc/MuonIdHelperSvc"};

◆ m_regionSelector

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

Definition at line 60 of file ClusterRoadDefiner.h.

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

◆ m_rWidth_RPC_Failed

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

Definition at line 57 of file ClusterRoadDefiner.h.

57{0};

◆ m_use_rpc

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

Definition at line 58 of file ClusterRoadDefiner.h.

58{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: