ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::InDetVertexSplitter Class Reference

#include <InDetVertexSplitter.h>

Inheritance diagram for InDet::InDetVertexSplitter:
Collaboration diagram for InDet::InDetVertexSplitter:

Public Member Functions

 InDetVertexSplitter (const std::string &name, ISvcLocator *pSvcLocator)
 Author: Peter V.
 ~InDetVertexSplitter ()
 Destructor - check up memory allocation delete any memory allocation on the heap.
StatusCode initialize ()
 Initialize initialize StoreGate.
StatusCode finalize ()
 Finalize - delete any memory allocation from the heap.
StatusCode execute ()
 Execute - on event by event.
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
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 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

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

Private Attributes

std::vector< std::string > m_trackKeys
bool m_isMatchedOdd
bool m_isUnmatchOdd
int m_addToVxMatched
int m_addToVxUnmatch
int m_eventN
std::string m_vertexContainerName
 containers to retrieve
std::string m_tpbContainerName
std::string m_trackContainerName
int m_maxVtx
bool m_priOnly
bool m_savetpb
DataObjIDColl m_extendedExtraObjects
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 31 of file InDetVertexSplitter.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ InDetVertexSplitter()

InDet::InDetVertexSplitter::InDetVertexSplitter ( const std::string & name,
ISvcLocator * pSvcLocator )

Author: Peter V.

Loscutoff (plosc.nosp@m.utof.nosp@m.f@gma.nosp@m.il.c.nosp@m.om February 2009 head file Constructor

switches to control the analysis through job options

Definition at line 36 of file InDetVertexSplitter.cxx.

37 :
38 AthAlgorithm(name, pSvcLocator),
39 m_isMatchedOdd(false),
40 m_isUnmatchOdd(false),
43 m_eventN(0){
44
46 declareProperty("TPBContainerName", m_tpbContainerName = "TrackParticleCandidate");
47 declareProperty("TrackContainerName", m_trackContainerName = "Tracks");
48 declareProperty("VertexContainerName", m_vertexContainerName = "VxPrimaryCandidate");
49 declareProperty("MaxVertexNumber", m_maxVtx = 1); //this should not be changed until a more robust handling of multi vertices is implemented
50 declareProperty("PrimaryOnly",m_priOnly = true); //this should not be changed presently
51 declareProperty("UseTrackParticleBase",m_savetpb = true); //is this needed?
52}
AthAlgorithm()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
std::string m_vertexContainerName
containers to retrieve

◆ ~InDetVertexSplitter()

InDet::InDetVertexSplitter::~InDetVertexSplitter ( )
default

Destructor - check up memory allocation delete any memory allocation on the heap.

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::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 }

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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< Algorithm > >::evtStore ( )
inlineinherited

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

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode InDet::InDetVertexSplitter::execute ( )

Execute - on event by event.

Definition at line 103 of file InDetVertexSplitter.cxx.

103 {
104
105 ATH_MSG_DEBUG("in execute()");
106
107 StatusCode sc = StatusCode::SUCCESS;
108
109 sc = split_vertices();
110 if (sc.isFailure()) {
111 ATH_MSG_ERROR("InDetVertexSplitter Failed");
112 return sc;
113 }
114
115 return sc;
116}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
static Double_t sc
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::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

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ finalize()

StatusCode InDet::InDetVertexSplitter::finalize ( )

Finalize - delete any memory allocation from the heap.

Definition at line 93 of file InDetVertexSplitter.cxx.

93 {
94 ATH_MSG_DEBUG("in finalize()");
95
96 return StatusCode::SUCCESS;
97
98}

◆ initialize()

StatusCode InDet::InDetVertexSplitter::initialize ( )

Initialize initialize StoreGate.

Definition at line 64 of file InDetVertexSplitter.cxx.

64 {
65
66 m_isMatchedOdd = false;
67 m_isUnmatchOdd = false;
70 m_eventN = 0;
71 std::stringstream ss;
72 for (int i = 1; i <=m_maxVtx; i++){
73 ss.str("");
74 ss << "odd_" << i << "_Tracks";
75 m_trackKeys.push_back(ss.str());
76 ss.str("");
77 ss << "even_" << i << "_Tracks";
78 m_trackKeys.push_back(ss.str());
79 ss.str("");
80 ss << "all_" << i << "_Tracks";
81 m_trackKeys.push_back(ss.str());
82 ss.str("");
83 }
84
85 ATH_MSG_INFO ("Initializing InDetVertexSplitter");
86
87 return StatusCode::SUCCESS;
88}
#define ATH_MSG_INFO(x)
static Double_t ss
std::vector< std::string > m_trackKeys

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

bool AthCommonMsg< Algorithm >::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< Algorithm > >::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< Algorithm > >::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< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ split_vertices()

StatusCode InDet::InDetVertexSplitter::split_vertices ( )
private

Definition at line 120 of file InDetVertexSplitter.cxx.

120 {
121
122 ATH_MSG_DEBUG("in split_vertices()");
123
124 StatusCode sc = StatusCode::SUCCESS;
125
126 const VxContainer* vtxTES=nullptr;
127 sc=evtStore()->retrieve( vtxTES, m_vertexContainerName);
128 if( sc.isFailure() || !vtxTES ) {
129 ATH_MSG_WARNING("No VxContainer container found in TDS tried " << m_vertexContainerName);
130 return StatusCode::SUCCESS;
131 }
132
133 ATH_MSG_DEBUG("Vertex successfully retrieved");
134
135 const Rec::TrackParticleContainer* tpbTES=nullptr;
136 const TrackCollection* trkTES=nullptr;
137
138 if (m_savetpb){
139 sc=evtStore()->retrieve( tpbTES, m_tpbContainerName);
140 if( sc.isFailure() || !tpbTES ) {
141 ATH_MSG_WARNING("No TrackParticleBase container found in TDS tried " << m_tpbContainerName);
142 return StatusCode::SUCCESS;
143 }
144 ATH_MSG_DEBUG("TrackParticleCandidate Collection successfully retrieved");
145 }
146 else {
147 sc=evtStore()->retrieve( trkTES, m_trackContainerName);
148 if( sc.isFailure() || !trkTES ) {
149 ATH_MSG_WARNING("No TrackCollection container found in TDS tried " << m_trackContainerName);
150 return StatusCode::SUCCESS;
151 }
152 ATH_MSG_DEBUG("TrackParticleCandidate Collection successfully retrieved");
153 }
154
155 std::map<std::string,TrackCollection*> trackmap;
156 std::map<std::string,Trk::TrackParticleBaseCollection*> tpbmap;
157
158// We need to create every container for each event, even if we don't write to them
159
160 for (const auto & thisKey : m_trackKeys){
161 TrackCollection* tempTracks = nullptr;
162 trackmap[thisKey] = tempTracks;
163 if (evtStore()->contains<TrackCollection>(thisKey) && (evtStore()->retrieve(trackmap[thisKey],thisKey)).isSuccess()){
164 //nop
165 } else {
166 trackmap[thisKey] = new TrackCollection;
167 }
168 }
169
170 for (const auto & thisKey : m_trackKeys){
171 Trk::TrackParticleBaseCollection* tempTpbs = nullptr;
172 tpbmap[thisKey] = tempTpbs;
174 (evtStore()->retrieve(tpbmap[thisKey],thisKey)).isSuccess()){
175 } else {
176 tpbmap[thisKey] = new Trk::TrackParticleBaseCollection;
177 }
178 }
179
180 if (m_savetpb and tpbTES){
181 //we loop over that list
182 for (const Rec::TrackParticle* tpb: *tpbTES){
183 const Trk::TrackParameters* trkPerigee = &(tpb->definingParameters());
184 bool trackmatched = false;
185 //we compare it to the tracks already associated with vertices
186
187 ATH_MSG_DEBUG("Found "<<vtxTES->size()<<" vertices");
188
189 int i_vtx = 0;
190 // We're relying here on the implicit sorting of vertices by sqrt(N_tracks)*Sum Pt_track^2
191 // This should pick out the most interesting N vertices
192 // Hopefully we have only 1 primary vertex, but if there is > 1 we can grab all of those too
193 std::stringstream sss;
194 std::string oeNameString;
195 oeNameString.reserve(20);
196 for (const Trk::VxCandidate* vtx : *vtxTES){
197 if ( (!m_priOnly || vtx->vertexType() == 1) && (i_vtx < m_maxVtx) ){
198 i_vtx++;
199 const std::vector<Trk::VxTrackAtVertex*> & vertexTracks = *vtx->vxTrackAtVertex();
200 ATH_MSG_DEBUG("parent vertex has "<<vertexTracks.size()<<" tracks, at position: "<<vtx->recVertex().position().x());
201 std::vector<Trk::VxTrackAtVertex*>::const_iterator tavI = vertexTracks.begin();
202 std::vector<Trk::VxTrackAtVertex*>::const_iterator tavIe= vertexTracks.end();
203 for (; tavI != tavIe; ++tavI){
204 const Trk::TrackParameters* vxTrkPerigee = (*tavI)->initialPerigee();
205 if (trkPerigee == vxTrkPerigee) {trackmatched = true;}
206 }
207 Trk::TrackParticleBase *trkCopy1 = new Trk::TrackParticleBase((*tpb));
208 Trk::TrackParticleBase *trkCopy2 = new Trk::TrackParticleBase((*tpb));
209 if (!trackmatched){
210 oeNameString.clear();
211 if (m_isUnmatchOdd) oeNameString = "odd";
212 if (!m_isUnmatchOdd) oeNameString = "even";
213 sss.str("");
214 sss << oeNameString << "_" << m_addToVxUnmatch << "_Tracks";
215 std::string oecontainerName = sss.str();
216 std::string allNameString = "all";
217 sss.str("");
218 sss << allNameString << "_" << m_addToVxUnmatch << "_Tracks";
219 std::string allcontainerName = sss.str();
220 ATH_MSG_DEBUG("found an unmatched trackparticlebase, giving it the key: "<< oecontainerName);
221 tpbmap[oecontainerName]->push_back(trkCopy1);
222 ATH_MSG_DEBUG("found an unmatched trackparticlebase, giving it the key: "<< allcontainerName);
223 tpbmap[allcontainerName]->push_back(trkCopy2);
227 }
228 if (trackmatched){
229 oeNameString.clear();
230 if (m_isMatchedOdd) oeNameString = "odd";
231 if (!m_isMatchedOdd) oeNameString = "even";
232 sss.str("");
233 sss << oeNameString << "_" << m_addToVxMatched << "_Tracks";
234 std::string oecontainerName = sss.str();
235 std::string allNameString = "all";
236 sss.str("");
237 sss << allNameString << "_" << m_addToVxMatched << "_Tracks";
238 std::string allcontainerName = sss.str();
239 ATH_MSG_DEBUG("found a matched trackparticlebase, giving it the key: "<< oecontainerName);
240 tpbmap[oecontainerName]->push_back(trkCopy1);
241 ATH_MSG_DEBUG("found a matched trackparticlebase, giving it the key: "<< allcontainerName);
242 tpbmap[allcontainerName]->push_back(trkCopy2);
246 }
247 }
248 }
249 }
250 }
251 if (!m_savetpb and trkTES){
252 TrackCollection::const_iterator trkItr = trkTES->begin();
253 TrackCollection::const_iterator trkItrE = trkTES->end();
254 //we loop over that list
255 for (; trkItr != trkItrE; ++trkItr){
256 const Trk::Perigee* trkPerigee = (*trkItr)->perigeeParameters();
257 bool trackmatched = false;
258 //we compare it to the tracks already associated with vertices
259
260 VxContainer::const_iterator vtxItr = vtxTES->begin();
261 VxContainer::const_iterator vtxItrE = vtxTES->end();
262
263 ATH_MSG_DEBUG("Found "<<vtxTES->size()<<" vertices");
264
265 int i_vtx = 0;
266 // We're relying here on the implicit sorting of vertices by sqrt(N_tracks)*Sum Pt_track^2
267 // This should pick out the most interesting N vertices
268 // Hopefully we have only 1 primary vertex, but if there is > 1 we can grab all of those too
269 std::stringstream sss;
270 std::string oeNameString;
271 oeNameString.reserve(20);
272 for (; vtxItr != vtxItrE; ++vtxItr){
273 if ( (!m_priOnly || (*vtxItr)->vertexType() == 1) && (i_vtx < m_maxVtx) ){
274 i_vtx++;
275 const std::vector<Trk::VxTrackAtVertex*> & vertexTracks = (*(*vtxItr)->vxTrackAtVertex());
276 ATH_MSG_DEBUG("parent vertex has "<<vertexTracks.size()<<" tracks, at position: "<<(*vtxItr)->recVertex().position().x());
277 std::vector<Trk::VxTrackAtVertex*>::const_iterator tavI = vertexTracks.begin();
278 std::vector<Trk::VxTrackAtVertex*>::const_iterator tavIe= vertexTracks.end();
279 for (; tavI != tavIe; ++tavI){
280 const Trk::TrackParameters* vxTrkPerigee = (*tavI)->initialPerigee();
281 if (trkPerigee == vxTrkPerigee) {trackmatched = true;}
282 }
283 Trk::Track *trkCopy1 = new Trk::Track((*(*trkItr)));
284 Trk::Track *trkCopy2 = new Trk::Track((*(*trkItr)));
285 if (!trackmatched){
286 oeNameString.clear();
287 if (m_isUnmatchOdd) oeNameString = "odd";
288 if (!m_isUnmatchOdd) oeNameString = "even";
289 sss.str("");
290 sss << oeNameString << "_" << m_addToVxUnmatch << "_Tracks";
291 std::string oecontainerName = sss.str();
292 std::string allNameString = "all";
293 sss.str("");
294 sss << allNameString << "_" << m_addToVxUnmatch << "_Tracks";
295 std::string allcontainerName = sss.str();
296 ATH_MSG_DEBUG("found an unmatched track, giving it the key: "<< oecontainerName);
297 trackmap[oecontainerName]->push_back(trkCopy1);
298 ATH_MSG_DEBUG("found an unmatched track, giving it the key: "<< allcontainerName);
299 trackmap[allcontainerName]->push_back(trkCopy2);
303 }
304 if (trackmatched){
305 oeNameString.clear();
306 if (m_isMatchedOdd) oeNameString = "odd";
307 if (!m_isMatchedOdd) oeNameString = "even";
308 sss.str("");
309 sss << oeNameString << "_" << m_addToVxMatched << "_Tracks";
310 std::string oecontainerName = sss.str();
311 std::string allNameString = "all";
312 sss.str("");
313 sss << allNameString << "_" << m_addToVxMatched << "_Tracks";
314 std::string allcontainerName = sss.str();
315 ATH_MSG_DEBUG("found a matched track, giving it the key: "<< oecontainerName);
316 trackmap[oecontainerName]->push_back(trkCopy1);
317 ATH_MSG_DEBUG("found a matched track, giving it the key: "<< allcontainerName);
318 trackmap[allcontainerName]->push_back(trkCopy2);
322 }
323 }
324 }
325 }
326 }
327 if (m_savetpb){
328 for (const auto & key: m_trackKeys){
329 if(evtStore()->record(tpbmap[key],key,false).isFailure() ){
330 ATH_MSG_ERROR("Could not save the "<< key);
331 }
332 }
333 } else {
334 for (const auto & key: m_trackKeys){
335 if(evtStore()->record(trackmap[key],key,false).isFailure() ){
336 ATH_MSG_ERROR("Could not save the "<< key);
337 }
338 }
339 }
340 ATH_MSG_DEBUG("split_vertices() succeeded");
341 m_eventN++;
342 return StatusCode::SUCCESS;
343}
#define ATH_MSG_WARNING(x)
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
Athena::TPCnvVers::Old VxContainer
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
size_type size() const noexcept
Returns the number of elements in the collection.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
DataVector< TrackParticleBase > TrackParticleBaseCollection
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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_addToVxMatched

int InDet::InDetVertexSplitter::m_addToVxMatched
private

Definition at line 48 of file InDetVertexSplitter.h.

◆ m_addToVxUnmatch

int InDet::InDetVertexSplitter::m_addToVxUnmatch
private

Definition at line 49 of file InDetVertexSplitter.h.

◆ m_detStore

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

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_eventN

int InDet::InDetVertexSplitter::m_eventN
private

Definition at line 50 of file InDetVertexSplitter.h.

◆ m_evtStore

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

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_isMatchedOdd

bool InDet::InDetVertexSplitter::m_isMatchedOdd
private

Definition at line 46 of file InDetVertexSplitter.h.

◆ m_isUnmatchOdd

bool InDet::InDetVertexSplitter::m_isUnmatchOdd
private

Definition at line 47 of file InDetVertexSplitter.h.

◆ m_maxVtx

int InDet::InDetVertexSplitter::m_maxVtx
private

Definition at line 58 of file InDetVertexSplitter.h.

◆ m_priOnly

bool InDet::InDetVertexSplitter::m_priOnly
private

Definition at line 59 of file InDetVertexSplitter.h.

◆ m_savetpb

bool InDet::InDetVertexSplitter::m_savetpb
private

Definition at line 60 of file InDetVertexSplitter.h.

◆ m_tpbContainerName

std::string InDet::InDetVertexSplitter::m_tpbContainerName
private

Definition at line 56 of file InDetVertexSplitter.h.

◆ m_trackContainerName

std::string InDet::InDetVertexSplitter::m_trackContainerName
private

Definition at line 57 of file InDetVertexSplitter.h.

◆ m_trackKeys

std::vector<std::string> InDet::InDetVertexSplitter::m_trackKeys
private

Definition at line 44 of file InDetVertexSplitter.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexContainerName

std::string InDet::InDetVertexSplitter::m_vertexContainerName
private

containers to retrieve

Definition at line 55 of file InDetVertexSplitter.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.


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