ATLAS Offline Software
Loading...
Searching...
No Matches
CaloClusterVertexFractionMaker Class Reference

#include <CaloClusterVertexFractionMaker.h>

Inheritance diagram for CaloClusterVertexFractionMaker:
Collaboration diagram for CaloClusterVertexFractionMaker:

Public Member Functions

 CaloClusterVertexFractionMaker (const std::string &type, const std::string &name, const IInterface *parent)
 ~CaloClusterVertexFractionMaker ()
virtual StatusCode initialize () override
virtual StatusCode execute (const EventContext &ctx, xAOD::CaloClusterContainer *theClusColl) const override
 Execute on an entire collection of clusters.
virtual StatusCode execute (xAOD::CaloClusterContainer *collection) final
 Execute on an entire collection of clusters.
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
 DeclareInterfaceID (CaloClusterCollectionProcessor, 1, 0)

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>

Static Private Member Functions

static double calculateDPhi (double phi1, double phi2)

Private Attributes

const double m_CALO_INNER_R
const double m_CALO_INNER_Z
double m_dRMatchMax
double m_dR2MatchMax
double m_maxClusterEta
ToolHandle< Trk::IExtrapolatorm_extrapolator
std::string m_vxContainerName
std::unique_ptr< Trk::Surfacem_cylinderSurface_atCaloEntrance
std::unique_ptr< Trk::Surfacem_discSurface_atCaloEntrance_positiveZ
std::unique_ptr< Trk::Surfacem_discSurface_atCaloEntrance_negativeZ
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 16 of file CaloClusterVertexFractionMaker.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

◆ CaloClusterVertexFractionMaker()

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

Definition at line 27 of file CaloClusterVertexFractionMaker.cxx.

28 : AthAlgTool(type, name, parent),
29 m_CALO_INNER_R(1450.), // millimeter
30 m_CALO_INNER_Z(3641.), // millimeter
31 m_dRMatchMax(0.1),
32 m_dR2MatchMax(0.), // will be overwritten in initialize
33 m_maxClusterEta(2.5),
34 m_extrapolator("Trk::Extrapolator/AtlasExtrapolator"),
35 m_vxContainerName("VxPrimaryCandidate")
36{
37 declareProperty ( "Extrapolator", m_extrapolator );
38 declareProperty ( "VxContainerName", m_vxContainerName );
39 declareProperty ( "dRMatchMax", m_dRMatchMax );
40 declareProperty ( "maxClusterEta", m_maxClusterEta );
41}
AthAlgTool()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ToolHandle< Trk::IExtrapolator > m_extrapolator

◆ ~CaloClusterVertexFractionMaker()

CaloClusterVertexFractionMaker::~CaloClusterVertexFractionMaker ( )

Definition at line 43 of file CaloClusterVertexFractionMaker.cxx.

44{
45 // all delete statements are in finalize (otherwise genconf crashes .. no idea why)
46}

Member Function Documentation

◆ calculateDPhi()

double CaloClusterVertexFractionMaker::calculateDPhi ( double phi1,
double phi2 )
staticprivate

Definition at line 230 of file CaloClusterVertexFractionMaker.cxx.

231{
232 double dPhi = std::fabs(phi1 - phi2);
233 if (dPhi > M_PI) dPhi = 2.*M_PI - dPhi;
234// if (dPhi > 3.1415926535897931) dPhi = 6.2831853071795862 - dPhi;
235 return dPhi;
236}
#define M_PI
bool dPhi(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)

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

◆ DeclareInterfaceID()

CaloClusterCollectionProcessor::DeclareInterfaceID ( CaloClusterCollectionProcessor ,
1 ,
0  )
inherited

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

◆ execute() [1/2]

StatusCode CaloClusterVertexFractionMaker::execute ( const EventContext & ctx,
xAOD::CaloClusterContainer * collection ) const
overridevirtual

Execute on an entire collection of clusters.

Parameters
collectionThe container of clusters. param ctx The event context.

usual complicated procedure to get from the vertex to the track(particle)

Implements CaloClusterCollectionProcessor.

Definition at line 73 of file CaloClusterVertexFractionMaker.cxx.

75{
76 const VxContainer* primcontainer(nullptr);
78 {
79 if ( evtStore()->retrieve ( primcontainer, m_vxContainerName ).isFailure() )
80 {
81 ATH_MSG_WARNING( "Could not retrieve collection " << m_vxContainerName << " in StoreGate, but contains<> says it is there." );
82 return StatusCode::FAILURE;
83 }
84 } else {
85 ATH_MSG_WARNING( "No collection " << m_vxContainerName << " in StoreGate." );
86 return StatusCode::FAILURE;
87 }
88
89 // loop over vertices, extrapolate tracks to calo, remember num tracks per vertex (for cluster vertex fraction calculation later)
90 std::vector<unsigned int> numTracksPerVertex(primcontainer->size()-1, 0);
91 std::vector<float> trkParticlePt_atOrigin;
92 std::vector<float> trkParticleEta_atCaloEntrance;
93 std::vector<float> trkParticlePhi_atCaloEntrance;
94
95 for (unsigned int v = 0 ; v < primcontainer->size()-1; ++v)
96 {
97 const std::vector<Trk::VxTrackAtVertex*>* vxTrackAtVertex = primcontainer->at(v)->vxTrackAtVertex();
98 for (std::vector<Trk::VxTrackAtVertex*>::const_iterator vxTrkItr = vxTrackAtVertex->begin(); vxTrkItr != vxTrackAtVertex->end(); ++vxTrkItr)
99 {
101 Trk::ITrackLink* trklink = (*vxTrkItr)->trackOrParticleLink();
102 Trk::LinkToTrackParticleBase* linkToTrackParticle = dynamic_cast<Trk::LinkToTrackParticleBase*>(trklink);
103 if (linkToTrackParticle != nullptr && linkToTrackParticle->isValid()) {
104 const Rec::TrackParticle* theTrackParticle = dynamic_cast<const Rec::TrackParticle*>(linkToTrackParticle->cachedElement());
105 if (theTrackParticle != nullptr)
106 {
107 const Trk::TrackParameters* trackParameters_atCaloEntrance(nullptr);
108 numTracksPerVertex.at(v)++;
109 trkParticlePt_atOrigin.push_back (theTrackParticle->pt()/1.e3);
110
111 const Trk::TrackParameters* lastTrackParametersInID(nullptr);
112 const std::vector< const Trk::TrackParameters * >& trackParametersVector = theTrackParticle->trackParameters();
113 if (trackParametersVector.size() > 1) lastTrackParametersInID = trackParametersVector.at(trackParametersVector.size()-2);
114 else lastTrackParametersInID = trackParametersVector.at(0); // the perigee
115
116 // very simplified way for now to decide to extrapolate a track to barrel or endcap
117 if (theTrackParticle->eta() >
118 1.35) // track most likely in endcap A, extrapolate track to a
119 // disc at ID exist in positive z direction
120 {
121 trackParameters_atCaloEntrance =
122 dynamic_cast<const Trk::AtaDisc*>(m_extrapolator->extrapolate(
123 ctx,
124 *lastTrackParametersInID,
127 true,
128 Trk::pion).release());
129 } else if (theTrackParticle->eta() <
130 -1.35) // track most likely in endcap C, extrapolate track
131 // to a disc at ID exist in negative z direction
132 {
133 trackParameters_atCaloEntrance =
134 dynamic_cast<const Trk::AtaDisc*>(m_extrapolator->extrapolate(
135 ctx,
136 *lastTrackParametersInID,
139 true,
140 Trk::pion).release());
141 } else // track is in barrel, extrapolate to cylinder at ID exit
142 {
143 trackParameters_atCaloEntrance =
144 dynamic_cast<const Trk::AtaCylinder*>(
145 m_extrapolator->extrapolate(ctx,
146 *lastTrackParametersInID,
149 true,
150 Trk::pion).release());
151 }
152 if (trackParameters_atCaloEntrance != nullptr) {
153 trkParticleEta_atCaloEntrance.push_back(trackParameters_atCaloEntrance->position().eta());
154 trkParticlePhi_atCaloEntrance.push_back(trackParameters_atCaloEntrance->position().phi());
155 ATH_MSG_DEBUG( "At calo entrance R(1150mm) " << *trackParameters_atCaloEntrance );
156 ATH_MSG_DEBUG( "TrkParticle eta/phi/pt[GeV] at calo " << trackParameters_atCaloEntrance->position().eta() << "\t"
157 << trackParameters_atCaloEntrance->position().phi() << "\t"
158 << trackParameters_atCaloEntrance->position().perp()/1.e3 );
159 delete trackParameters_atCaloEntrance;
160 } else {
161 trkParticleEta_atCaloEntrance.push_back(999.);
162 trkParticlePhi_atCaloEntrance.push_back(999.);
163 }
164 }
165 }
166 }
167 }
168
169 double fabs_Dphi(9999.);
170 double fabs_Deta(9999.);
171 double dR2(9999.);
172 std::vector<float> sumPtOfMatchedTracksPerVertex(numTracksPerVertex.size(), 0.);
173
174 for (xAOD::CaloClusterContainer::iterator clItr = caloClusterContainer->begin(); clItr != caloClusterContainer->end(); ++clItr)
175 {
176 xAOD::CaloCluster* theCluster = (*clItr);
177// double cl_theta = 2.*atan(exp(-theCluster->eta()));
178// double cl_et = theCluster->e()/1.e3*sin(cl_theta);
179 if (std::fabs(theCluster->eta()) < 2.5)
180 {
181 // Calculation of Cluster Vertex Fraction
182 std::vector<unsigned int>::iterator numTrksPerVertexItr = numTracksPerVertex.begin();
183 std::vector<unsigned int>::iterator numTrksPerVertexItrE = numTracksPerVertex.end();
184 unsigned int vertexCounter(0);
185 unsigned int totalTrackCounter(0);
186 for ( ; numTrksPerVertexItr != numTrksPerVertexItrE ; ++numTrksPerVertexItr, vertexCounter++)
187 {
188 sumPtOfMatchedTracksPerVertex.at(vertexCounter) = 0.;
189 for (unsigned int track = 0; track < (*numTrksPerVertexItr); ++track, totalTrackCounter++)
190 {
191 if (trkParticleEta_atCaloEntrance.at(totalTrackCounter) < 900.) // no need to check phi as well (not extrapolated eta is set to 999.)
192 {
193 fabs_Dphi = calculateDPhi(trkParticlePhi_atCaloEntrance.at(totalTrackCounter), theCluster->phi());
194 fabs_Deta = fabs(trkParticleEta_atCaloEntrance.at(totalTrackCounter) - theCluster->eta());
195 dR2 = fabs_Deta * fabs_Deta + fabs_Dphi * fabs_Dphi;
196 if (dR2 < m_dR2MatchMax) {
197 sumPtOfMatchedTracksPerVertex.at(vertexCounter) += trkParticlePt_atOrigin.at(totalTrackCounter);
198 }
199 }
200 }
201 }
202 double totalSumPtOfMatchedTracksWhichWereAlsoUsedInAVertex(0.);
203 for (unsigned int mTpV = 0; mTpV < sumPtOfMatchedTracksPerVertex.size(); ++mTpV) {
204 totalSumPtOfMatchedTracksWhichWereAlsoUsedInAVertex += sumPtOfMatchedTracksPerVertex.at(mTpV);
205 //std::cout << "CaloClusterVertexFractionMakerAOD: " << mTpV << "\t" << sumPtOfMatchedTracksPerVertex.at(mTpV) << std::endl;
206 }
207
208 double cvf(-1.);
209 double ncvf(-1.);
210 if (totalSumPtOfMatchedTracksWhichWereAlsoUsedInAVertex > 0.)
211 {
212 double sumPtInPrimary = sumPtOfMatchedTracksPerVertex.at(0);
213 cvf = sumPtOfMatchedTracksPerVertex.at(0)/totalSumPtOfMatchedTracksWhichWereAlsoUsedInAVertex;
214 std::sort(sumPtOfMatchedTracksPerVertex.begin(), sumPtOfMatchedTracksPerVertex.end());
215 double highestSumPt = sumPtOfMatchedTracksPerVertex.at(sumPtOfMatchedTracksPerVertex.size()-1);
216 if (highestSumPt > 0) ncvf = sumPtInPrimary/highestSumPt;
217 }
218 //std::cout << "CaloClusterVertexFractionMakerAOD: cvf " << cvf << "\tncvf " << ncvf << std::endl;
219 theCluster->insertMoment(xAOD::CaloCluster::VERTEX_FRACTION, cvf); // false at the end writes the moment into the moment store of the cluster directly (the moment is then written to AOD as well)
220 theCluster->insertMoment(xAOD::CaloCluster::NVERTEX_FRACTION, ncvf); // false at the end writes the moment into the moment store of the cluster directly (the moment is then written to AOD as well)
221 } else {
222 theCluster->insertMoment(xAOD::CaloCluster::VERTEX_FRACTION, -1.); // false at the end writes the moment into the moment store of the cluster directly (the moment is then written to AOD as well)
223 theCluster->insertMoment(xAOD::CaloCluster::NVERTEX_FRACTION, -1.); // false at the end writes the moment into the moment store of the cluster directly (the moment is then written to AOD as well)
224 }
225 }
226
227 return StatusCode::SUCCESS;
228}
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Old VxContainer
ServiceHandle< StoreGateSvc > & evtStore()
std::unique_ptr< Trk::Surface > m_discSurface_atCaloEntrance_positiveZ
std::unique_ptr< Trk::Surface > m_cylinderSurface_atCaloEntrance
std::unique_ptr< Trk::Surface > m_discSurface_atCaloEntrance_negativeZ
static double calculateDPhi(double phi1, double phi2)
DataModel_detail::iterator< DataVector > iterator
Definition DataVector.h:842
virtual double pt() const =0
transverse momentum
virtual double eta() const =0
pseudo rapidity
const std::vector< const TrackParameters * > & trackParameters() const
Returns the track parameters.
virtual double eta() const
The pseudorapidity ( ) of the particle.
void insertMoment(MomentType type, double value)
virtual double phi() const
The azimuthal angle ( ) of the particle.
@ VERTEX_FRACTION
Vertex fraction of this cluster wrt.
@ NVERTEX_FRACTION
slightly updated vertex fraction more pile up independent (similar to nJVF)
static std::string release
Definition computils.h:50
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
ParametersT< TrackParametersDim, Charged, DiscSurface > AtaDisc
@ alongMomentum
ParametersBase< TrackParametersDim, Charged > TrackParameters
ParametersT< TrackParametersDim, Charged, CylinderSurface > AtaCylinder
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.

◆ execute() [2/2]

virtual StatusCode CaloClusterCollectionProcessor::execute ( xAOD::CaloClusterContainer * collection)
inlinefinalvirtual

Execute on an entire collection of clusters.

Parameters
collectionThe container of clusters. (deprecated)

Reimplemented from CaloClusterCollectionProcessor.

Definition at line 50 of file CaloClusterCollectionProcessor.h.

51 {
52 return execute (Gaudi::Hive::currentContext(), collection);
53 }
virtual StatusCode execute(const EventContext &ctx, xAOD::CaloClusterContainer *theClusColl) const override
Execute on an entire collection of clusters.

◆ 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 CaloClusterVertexFractionMaker::initialize ( )
overridevirtual

Definition at line 48 of file CaloClusterVertexFractionMaker.cxx.

49{
50 ATH_MSG_INFO( "Initializing " << name() );
51 ATH_MSG_INFO( "Only considering clusters with |eta| < " << m_maxClusterEta << "\tdRMatchMax Cluster-Track = " << m_dRMatchMax );
52
54
55 /* Get the extrapolator */
56 ATH_CHECK( m_extrapolator.retrieve() );
57 ATH_MSG_INFO( "Retrieved tool " << m_extrapolator );
58
59 // a surface at the entrance to the calorimeter
60// HepGeom::TranslateZ3D* translateAlongPositiveZ = new HepGeom::TranslateZ3D(m_CALO_INNER_Z);
61// HepGeom::TranslateZ3D* translateAlongNegativeZ = new HepGeom::TranslateZ3D(-m_CALO_INNER_Z);
62 Amg::Transform3D translateAlongPositiveZ = Amg::Transform3D(Amg::Vector3D(0.,0.,m_CALO_INNER_Z));
63 Amg::Transform3D translateAlongNegativeZ = Amg::Transform3D(Amg::Vector3D(0.,0.,m_CALO_INNER_Z));
64
65 m_cylinderSurface_atCaloEntrance = std::make_unique<Trk::CylinderSurface>(m_CALO_INNER_R, 8000.);
66 m_discSurface_atCaloEntrance_positiveZ = std::make_unique<Trk::DiscSurface>(translateAlongPositiveZ, 0., 10000.);
67 m_discSurface_atCaloEntrance_negativeZ = std::make_unique<Trk::DiscSurface>(translateAlongNegativeZ, 0., 10000.);
68
69 return StatusCode::SUCCESS;
70}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D

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

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

const double CaloClusterVertexFractionMaker::m_CALO_INNER_R
private

Definition at line 30 of file CaloClusterVertexFractionMaker.h.

◆ m_CALO_INNER_Z

const double CaloClusterVertexFractionMaker::m_CALO_INNER_Z
private

Definition at line 31 of file CaloClusterVertexFractionMaker.h.

◆ m_cylinderSurface_atCaloEntrance

std::unique_ptr<Trk::Surface> CaloClusterVertexFractionMaker::m_cylinderSurface_atCaloEntrance
private

Definition at line 43 of file CaloClusterVertexFractionMaker.h.

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

std::unique_ptr<Trk::Surface> CaloClusterVertexFractionMaker::m_discSurface_atCaloEntrance_negativeZ
private

Definition at line 45 of file CaloClusterVertexFractionMaker.h.

◆ m_discSurface_atCaloEntrance_positiveZ

std::unique_ptr<Trk::Surface> CaloClusterVertexFractionMaker::m_discSurface_atCaloEntrance_positiveZ
private

Definition at line 44 of file CaloClusterVertexFractionMaker.h.

◆ m_dR2MatchMax

double CaloClusterVertexFractionMaker::m_dR2MatchMax
private

Definition at line 34 of file CaloClusterVertexFractionMaker.h.

◆ m_dRMatchMax

double CaloClusterVertexFractionMaker::m_dRMatchMax
private

Definition at line 33 of file CaloClusterVertexFractionMaker.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_extrapolator

ToolHandle< Trk::IExtrapolator > CaloClusterVertexFractionMaker::m_extrapolator
private

Definition at line 38 of file CaloClusterVertexFractionMaker.h.

◆ m_maxClusterEta

double CaloClusterVertexFractionMaker::m_maxClusterEta
private

Definition at line 35 of file CaloClusterVertexFractionMaker.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.

◆ m_vxContainerName

std::string CaloClusterVertexFractionMaker::m_vxContainerName
private

Definition at line 41 of file CaloClusterVertexFractionMaker.h.


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