ATLAS Offline Software
Loading...
Searching...
No Matches
Analysis::MSVVariablesFactory Class Reference

#include <MSVVariablesFactory.h>

Inheritance diagram for Analysis::MSVVariablesFactory:
Collaboration diagram for Analysis::MSVVariablesFactory:

Public Member Functions

 MSVVariablesFactory (const std::string &name, const std::string &n, const IInterface *p)
virtual ~MSVVariablesFactory ()=default
virtual StatusCode initialize () override
virtual StatusCode finalize () override
virtual StatusCode fillMSVVariables (const xAOD::Jet &, xAOD::BTagging *BTag, const Trk::VxSecVKalVertexInfo *myInfoVKal, xAOD::VertexContainer *btagVertex, const xAOD::Vertex &PV, std::string basename) const override
virtual StatusCode createMSVContainer (const xAOD::Jet &, const Trk::VxSecVKalVertexInfo *myInfoVKal, xAOD::VertexContainer *btagVertex, const xAOD::Vertex &PV) const override
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

Static Public Member Functions

static const InterfaceID & interfaceID ()

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

double get3DSignificance (const xAOD::Vertex *priVertex, std::vector< const xAOD::Vertex * > &secVertex, const Amg::Vector3D jetDirection) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

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 23 of file MSVVariablesFactory.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

◆ MSVVariablesFactory()

Analysis::MSVVariablesFactory::MSVVariablesFactory ( const std::string & name,
const std::string & n,
const IInterface * p )

Definition at line 39 of file MSVVariablesFactory.cxx.

41 :
42 AthAlgTool(name, n,p)
43 {
44 declareInterface<IMSVVariablesFactory>(this);
45 }
AthAlgTool()
Default constructor:

◆ ~MSVVariablesFactory()

virtual Analysis::MSVVariablesFactory::~MSVVariablesFactory ( )
virtualdefault

Member Function Documentation

◆ createMSVContainer()

StatusCode Analysis::MSVVariablesFactory::createMSVContainer ( const xAOD::Jet & myJet,
const Trk::VxSecVKalVertexInfo * myInfoVKal,
xAOD::VertexContainer * btagVertex,
const xAOD::Vertex & PV ) const
overridevirtual

Implements Analysis::IMSVVariablesFactory.

Definition at line 57 of file MSVVariablesFactory.cxx.

59 {
60
61 Amg::Vector3D jet_V3(myJet.p4().Px(), myJet.p4().Py(), myJet.p4().Pz());
62 float jetenergy=0.;
63 const xAOD::Vertex* priVtx = &PrimaryVtx;
64 std::vector< ElementLink< xAOD::VertexContainer > > MSVVertexLinks;
65 const std::vector<xAOD::Vertex*> myVertices = myVertexInfoVKal->vertices();
66 if(myVertices.empty()){
67 ATH_MSG_DEBUG("#BTAG# no MSV vertices...fill default values only... ");
69 VertexContainer->push_back(vertex);
77 return StatusCode::SUCCESS;
78 }
79
80 jetenergy = myVertexInfoVKal->energyTrkInJet();
81
82 for (const auto& vertex : myVertexInfoVKal->vertices()){
83 VertexContainer->push_back(vertex);
84 //additional info per vertex
85 double sumpx = 0.0;
86 double sumpy = 0.0;
87 double sumpz = 0.0;
88 double sume = 0.0;
89 const std::vector<ElementLink<xAOD::TrackParticleContainer> > myTrackLinks = vertex->trackParticleLinks();
90 if (myTrackLinks.empty()) {
91 ATH_MSG_WARNING("#BTAG# No Track Links attached to the track at the sec vertex... ");
92 }
93 int npsec = 0;
94 const std::vector<Trk::VxTrackAtVertex> myTracks=vertex->vxTrackAtVertex();
95 if (!myTracks.empty()) {
96 npsec=myTracks.size();
97 for (const auto& track : myTracks) {
98 const Trk::Perigee* perigee = dynamic_cast<const Trk::Perigee*>(track.perigeeAtVertex());
99 if(perigee){
100 sumpx += perigee->momentum().x();
101 sumpy += perigee->momentum().y();
102 sumpz += perigee->momentum().z();
103 sume += std::hypot(perigee->momentum().mag(), ParticleConstants::chargedPionMassInMeV);
104 }else{
105 ATH_MSG_WARNING("#BTAG# perigee for VxTrackAtVertex not found");
106 }
107 }
108 }
109
110 CLHEP::HepLorentzVector vtxp4(sumpx,sumpy,sumpz,sume);
111 float efrac = (jetenergy>0) ? vtxp4.e()/jetenergy : 0;
112 xAOD::SecVtxHelper::setVertexMass(vertex, vtxp4.m());
114 xAOD::SecVtxHelper::setVtxNtrk(vertex, npsec);
115 xAOD::SecVtxHelper::setVtxpt(vertex, vtxp4.perp());
116 xAOD::SecVtxHelper::setVtxeta(vertex, vtxp4.eta());
117 xAOD::SecVtxHelper::setVtxphi(vertex, vtxp4.phi());
118
119 ATH_MSG_DEBUG("#BTAG# mass per vertex = "<<vtxp4.m());
120 double localdistnrm = 0;
121 std::vector<const xAOD::Vertex*> vecVtxHolder;
122 vecVtxHolder.push_back(vertex);
123
124 if (priVtx) {
125 ATH_MSG_DEBUG("Factory PVX x = " << priVtx->x() << " y = " << priVtx->y() << " z = " << priVtx->z());
126 localdistnrm = get3DSignificance(priVtx, vecVtxHolder, jet_V3);
127 } else {
128 ATH_MSG_WARNING("#BTAG# Tagging requested, but no primary vertex supplied.");
129 }
130 xAOD::SecVtxHelper::setVtxnormDist(vertex, localdistnrm);
131 //track links,
132 vertex->setTrackParticleLinks(myTrackLinks);
133
134 } //end loop vertexcontainer
135
136 return StatusCode::SUCCESS;
137 }
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
double get3DSignificance(const xAOD::Vertex *priVertex, std::vector< const xAOD::Vertex * > &secVertex, const Amg::Vector3D jetDirection) const
const Amg::Vector3D & momentum() const
Access method for the momentum.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
float z() const
Returns the z position.
float y() const
Returns the y position.
float x() const
Returns the x position.
Eigen::Matrix< double, 3, 1 > Vector3D
constexpr double chargedPionMassInMeV
the mass of the charged pion (in MeV)
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
void setVtxeta(xAOD::Vertex *, float value)
void setVertexMass(xAOD::Vertex *, float value)
void setVtxNtrk(xAOD::Vertex *, int value)
void setVtxpt(xAOD::Vertex *, float value)
void setVtxnormDist(xAOD::Vertex *, float value)
void setEnergyFraction(xAOD::Vertex *, float value)
void setVtxphi(xAOD::Vertex *, float value)
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ 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>

◆ 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

◆ fillMSVVariables()

StatusCode Analysis::MSVVariablesFactory::fillMSVVariables ( const xAOD::Jet & myJet,
xAOD::BTagging * BTag,
const Trk::VxSecVKalVertexInfo * myInfoVKal,
xAOD::VertexContainer * btagVertex,
const xAOD::Vertex & PV,
std::string basename ) const
overridevirtual

Implements Analysis::IMSVVariablesFactory.

Definition at line 139 of file MSVVariablesFactory.cxx.

143 {
144
145 Amg::Vector3D jet_V3(myJet.p4().Px(), myJet.p4().Py(), myJet.p4().Pz());
146 int nvsec = 0;
147 float jetenergy = 0.;
148 int n2t = 0;
149 float distnrm = 0.;
150 const xAOD::Vertex* priVtx = &PrimaryVtx;
151 std::vector< ElementLink< xAOD::VertexContainer > > MSVVertexLinks;
152 const std::vector<xAOD::Vertex*> myVertices = myVertexInfoVKal->vertices();
153 if(myVertices.empty()){
154 ATH_MSG_DEBUG("#BTAG# no MSV vertices...fill default values only... ");
155 BTag->setVariable<int>(basename, "N2Tpair", n2t);
156 BTag->setVariable<float>(basename, "energyTrkInJet", jetenergy);
157 BTag->setVariable<int>(basename, "nvsec", nvsec);
158 BTag->setVariable<float>(basename, "normdist", distnrm);
159 BTag->setVariable<std::vector<ElementLink<xAOD::VertexContainer> > >(basename, "vertices", MSVVertexLinks);
160 BTag->setDynVxELName(basename, "vertices");
162 VertexContainer->push_back(vertex);
166 xAOD::SecVtxHelper::setVtxpt(vertex, -9.);
170 return StatusCode::SUCCESS;
171 }
172
173 jetenergy = myVertexInfoVKal->energyTrkInJet();
174 n2t = myVertexInfoVKal->n2trackvertices();
175 BTag->setVariable<int>(basename, "N2Tpair", n2t);
176 BTag->setVariable<float>(basename, "energyTrkInJet", jetenergy);
177
178 std::vector<const xAOD::Vertex*> vecVertices;
179 for (const auto& vertex : myVertexInfoVKal->vertices()) {
180 VertexContainer->push_back(vertex);
181 //additional info per vertex
182 vecVertices.push_back(vertex);
183 double sumpx = 0.0;
184 double sumpy = 0.0;
185 double sumpz = 0.0;
186 double sume = 0.0;
187 const std::vector<ElementLink<xAOD::TrackParticleContainer> > myTrackLinks = vertex->trackParticleLinks();
188 if (myTrackLinks.empty()) {
189 ATH_MSG_WARNING("#BTAG# No Track Links attached to the track at the sec vertex... ");
190 }
191 int npsec = 0;
192 const std::vector<Trk::VxTrackAtVertex> myTracks=vertex->vxTrackAtVertex();
193 if (!myTracks.empty()) {
194 npsec=myTracks.size();
195 for (const auto& track : myTracks) {
196 const Trk::Perigee* perigee = dynamic_cast<const Trk::Perigee*>(track.perigeeAtVertex());
197 if(perigee){
198 sumpx += perigee->momentum().x();
199 sumpy += perigee->momentum().y();
200 sumpz += perigee->momentum().z();
201 sume += std::hypot(perigee->momentum().mag(), ParticleConstants::chargedPionMassInMeV);
202 }else{
203 ATH_MSG_WARNING("#BTAG# perigee for VxTrackAtVertex not found");
204 }
205 }
206 }
207
208 CLHEP::HepLorentzVector vtxp4(sumpx,sumpy,sumpz,sume);
209 float efrac = (jetenergy>0) ? vtxp4.e()/jetenergy : 0;
210 xAOD::SecVtxHelper::setVertexMass(vertex, vtxp4.m());
212 xAOD::SecVtxHelper::setVtxNtrk(vertex, npsec);
213 xAOD::SecVtxHelper::setVtxpt(vertex, vtxp4.perp());
214 xAOD::SecVtxHelper::setVtxeta(vertex, vtxp4.eta());
215 xAOD::SecVtxHelper::setVtxphi(vertex, vtxp4.phi());
216
217 ATH_MSG_DEBUG("#BTAG# mass per vertex = "<<vtxp4.m());
218 double localdistnrm = 0;
219 std::vector<const xAOD::Vertex*> vecVtxHolder;
220 vecVtxHolder.push_back(vertex);
221
222
223 if (priVtx) {
224 ATH_MSG_DEBUG("Factory PVX x = " << priVtx->x() << " y = " << priVtx->y() << " z = " << priVtx->z());
225 localdistnrm = get3DSignificance(priVtx, vecVtxHolder, jet_V3);
226 } else {
227 ATH_MSG_WARNING("#BTAG# Tagging requested, but no primary vertex supplied.");
228 }
229 xAOD::SecVtxHelper::setVtxnormDist(vertex, localdistnrm);
230 //track links,
231 vertex->setTrackParticleLinks(myTrackLinks);
232
233 ElementLink< xAOD::VertexContainer> linkBTagVertex;
234 linkBTagVertex.toContainedElement(*VertexContainer, vertex);
235 MSVVertexLinks.push_back(linkBTagVertex);
236 } //end loop vertexcontainer
237
238 BTag->setVariable<std::vector<ElementLink<xAOD::VertexContainer> > >(basename, "vertices", MSVVertexLinks);
239 BTag->setDynVxELName(basename, "vertices");
240
241 if (priVtx) {
242 distnrm = get3DSignificance(priVtx, vecVertices, jet_V3);
243 } else {
244 ATH_MSG_WARNING("#BTAG# Tagging requested, but no primary vertex supplied.");
245 distnrm=0.;
246 }
247 nvsec = vecVertices.size();
248 BTag->setVariable<int>(basename, "nvsec", nvsec);
249 BTag->setVariable<float>(basename, "normdist", distnrm);
250
251 return StatusCode::SUCCESS;
252
253 }
@ BTag
The object is a b-tagging object.
Definition ObjectType.h:60
std::string basename(std::string name)
Definition utils.cxx:207

◆ finalize()

StatusCode Analysis::MSVVariablesFactory::finalize ( )
overridevirtual

Implements Analysis::IMSVVariablesFactory.

Definition at line 52 of file MSVVariablesFactory.cxx.

52 {
53 ATH_MSG_DEBUG(" Finalization of MSVVariablesFactory succesfull");
54 return StatusCode::SUCCESS;
55 }

◆ get3DSignificance()

double Analysis::MSVVariablesFactory::get3DSignificance ( const xAOD::Vertex * priVertex,
std::vector< const xAOD::Vertex * > & secVertex,
const Amg::Vector3D jetDirection ) const
private

Definition at line 255 of file MSVVariablesFactory.cxx.

258 {
259
260 if(!secVertex.size()) return 0;
261 std::vector<Amg::Vector3D> positions;
262 std::vector<AmgSymMatrix(3)> weightMatrices;
263 Amg::Vector3D weightTimesPosition(0.,0.,0.);
264 AmgSymMatrix(3) sumWeights;
265 sumWeights.setZero();
266
267 for (const auto& vertex : secVertex) {
268 positions.push_back(vertex->position());
269 weightMatrices.push_back(vertex->covariancePosition().inverse());
270 weightTimesPosition += weightMatrices.back() * positions.back();
271 sumWeights += weightMatrices.back();
272 }
273
274 bool invertible;
275 AmgSymMatrix(3) meanCovariance;
276 meanCovariance.setZero();
277 sumWeights.computeInverseWithCheck(meanCovariance, invertible);
278 if (!invertible) {
279 ATH_MSG_WARNING("#BTAG# Could not invert sum of sec vtx matrices");
280 return 0.;
281 }
282 Amg::Vector3D meanPosition = meanCovariance * weightTimesPosition;
283 AmgSymMatrix(3) covariance = meanCovariance + priVertex->covariancePosition();
284
285 double Lx = meanPosition[0]-priVertex->position().x();
286 double Ly = meanPosition[1]-priVertex->position().y();
287 double Lz = meanPosition[2]-priVertex->position().z();
288
289 const double decaylength = sqrt(Lx*Lx + Ly*Ly + Lz*Lz);
290 const double inv_decaylength = 1. / decaylength;
291 double dLdLx = Lx * inv_decaylength;
292 double dLdLy = Ly * inv_decaylength;
293 double dLdLz = Lz * inv_decaylength;
294 double decaylength_err = sqrt(dLdLx*dLdLx*covariance(0,0) +
295 dLdLy*dLdLy*covariance(1,1) +
296 dLdLz*dLdLz*covariance(2,2) +
297 2.*dLdLx*dLdLy*covariance(0,1) +
298 2.*dLdLx*dLdLz*covariance(0,2) +
299 2.*dLdLy*dLdLz*covariance(1,2));
300
301 double decaylength_significance = 0.;
302 if (decaylength_err != 0.) decaylength_significance = decaylength/decaylength_err;
303 double L_proj_jetDir = jetDirection.x()*Lx + jetDirection.y()*Ly + jetDirection.z()*Lz;
304 if (L_proj_jetDir < 0.) decaylength_significance *= -1.;
305
306 return decaylength_significance;
307
308 }
#define AmgSymMatrix(dim)
if(febId1==febId2)
#define y
#define x
#define z

◆ initialize()

StatusCode Analysis::MSVVariablesFactory::initialize ( )
overridevirtual

Implements Analysis::IMSVVariablesFactory.

Definition at line 47 of file MSVVariablesFactory.cxx.

47 {
48 ATH_MSG_DEBUG(" Initialization of MSVVariablesFactory succesfull");
49 return StatusCode::SUCCESS;
50 }

◆ 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.

◆ interfaceID()

const InterfaceID & Analysis::IMSVVariablesFactory::interfaceID ( )
inlinestaticinherited

Definition at line 49 of file IMSVVariablesFactory.h.

49{ return IID_IMSVVariablesFactory; };
static const InterfaceID IID_IMSVVariablesFactory("Analysis::IMSVVariablesFactory", 1, 0)

◆ 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 }

◆ 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_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: