ATLAS Offline Software
Loading...
Searching...
No Matches
Trk::RecMomentumQualityValidation Class Reference

This is for automatising the 100%, 50%, 20%, 10% tables for (indet) tracks. More...

#include <RecMomentumQualityValidation.h>

Inheritance diagram for Trk::RecMomentumQualityValidation:
Collaboration diagram for Trk::RecMomentumQualityValidation:

Public Member Functions

 RecMomentumQualityValidation (const std::string &name, ISvcLocator *pSvcLocator)
 Standard Athena-Algorithm Constructor.
 ~RecMomentumQualityValidation ()
 Default Destructor.
StatusCode initialize ()
 standard Athena-Algorithm method
StatusCode execute ()
 standard Athena-Algorithm method
StatusCode finalize ()
 standard Athena-Algorithm method
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

enum  StatIndex { iAll = 0 , iBarrel = 1 , iTransi = 2 , iEndcap = 3 }
 define ranges in eta More...
typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

void printTable () const
 method: make the output table
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 void monitorTrackFits (std::vector< unsigned int > &, const double &)

Private Attributes

std::string m_inputTrackCollection
 properties from JobOptions:
std::string m_trackTruthCollection
 job option: the truth track collection name
ToolHandle< Trk::ITruthToTrackm_truthToTrack
 Tool handle to Trk::ITruthToTrack tool.
ToolHandle< Trk::ITrackSelectorToolm_trackSelector
 Tool handle to Trk::ITrackSelectorTool.
const AtlasDetectorIDm_idHelper
std::vector< unsigned int > m_nHundred
 counters
std::vector< unsigned int > m_nFifty
std::vector< unsigned int > m_nTwenty
std::vector< unsigned int > m_nTen
std::vector< unsigned int > m_nFakeOrLost
std::vector< unsigned int > m_tHundred
std::vector< unsigned int > m_tFifty
std::vector< unsigned int > m_tTwenty
std::vector< unsigned int > m_tTen
std::vector< unsigned int > m_tFakeOrLost
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

This is for automatising the 100%, 50%, 20%, 10% tables for (indet) tracks.

Author
Wolfgang Liebig (Wolfgang.Liebig -at- cern.ch)

Definition at line 32 of file RecMomentumQualityValidation.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Member Enumeration Documentation

◆ StatIndex

Constructor & Destructor Documentation

◆ RecMomentumQualityValidation()

Trk::RecMomentumQualityValidation::RecMomentumQualityValidation ( const std::string & name,
ISvcLocator * pSvcLocator )

Standard Athena-Algorithm Constructor.

Definition at line 22 of file RecMomentumQualityValidation.cxx.

22 :
23 AthAlgorithm(name,pSvcLocator),
24 m_inputTrackCollection("Tracks"),
25 m_trackTruthCollection("TrackTruthCollection"),
26 m_truthToTrack("Trk::TruthToTrack/TruthToTrack"),
27 m_trackSelector(""), //InDet::InDetTrackSelectorTool/InDetTrackSelectorTool),
28 m_idHelper(nullptr),
31{
32
33 // Get Parameter values from JobOptions file
34 declareProperty("InputTrackCollection", m_inputTrackCollection);
35 declareProperty("TrackTruthCollection", m_trackTruthCollection);
36 declareProperty("TruthToTrackTool", m_truthToTrack);
37 declareProperty("TrackSelectorTool", m_trackSelector);
38
39}
AthAlgorithm()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
std::string m_inputTrackCollection
properties from JobOptions:
std::vector< unsigned int > m_nHundred
counters
ToolHandle< Trk::ITrackSelectorTool > m_trackSelector
Tool handle to Trk::ITrackSelectorTool.
std::string m_trackTruthCollection
job option: the truth track collection name
ToolHandle< Trk::ITruthToTrack > m_truthToTrack
Tool handle to Trk::ITruthToTrack tool.

◆ ~RecMomentumQualityValidation()

Trk::RecMomentumQualityValidation::~RecMomentumQualityValidation ( )
default

Default Destructor.

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 Trk::RecMomentumQualityValidation::execute ( )

standard Athena-Algorithm method

Definition at line 93 of file RecMomentumQualityValidation.cxx.

94{
95
96 // Retrieving the Trackcollection specified via m_inputTrackCollection
97 StatusCode sc = StatusCode::SUCCESS;
98 const TrackCollection* trackCollection = nullptr;
99
100 if (!m_inputTrackCollection.empty()) {
101 sc = evtStore()->retrieve(trackCollection,m_inputTrackCollection);
102 if (sc.isFailure())
103 ATH_MSG_ERROR( "TrackCollection "<<m_inputTrackCollection<<" not found!" );
104 else
105 ATH_MSG_VERBOSE("TrackCollection " << m_inputTrackCollection<<" found." );
106
107 } else {
108 ATH_MSG_ERROR("No Track collection given!" );
109 return StatusCode::FAILURE;
110 }
111
112 const TrackTruthCollection* trackTruthCollection = nullptr;
113 if (!m_trackTruthCollection.empty()) {
114 sc = evtStore()->retrieve(trackTruthCollection, m_trackTruthCollection);
115 if (sc.isFailure())
116 ATH_MSG_ERROR( "TruthCollection "<<m_trackTruthCollection<<" not found!" );
117 else
118 ATH_MSG_VERBOSE("TruthCollection " << m_trackTruthCollection<<" found." );
119 } else {
120 ATH_MSG_ERROR("No Truth collection given!" );
121 return StatusCode::FAILURE;
122 }
123
124 // Looping over the tracks of retrieved trackcollection
125 TrackCollection::const_iterator trackIterator = trackCollection->begin();
126 for (;trackIterator!=trackCollection->end();++trackIterator) {
127
128 if (!m_trackSelector.empty() && (*trackIterator)!=nullptr &&
129 m_trackSelector->decision(**trackIterator) ) {
130 const Trk::TrackParameters* generatedTrackPerigee(nullptr);
131
132 // find matching truth particle
133 const TrackTruth* trackTruth = nullptr;
134 HepMC::ConstGenParticlePtr genParticle{nullptr};
135 TrackTruthCollection::const_iterator truthIterator
136 = trackTruthCollection->find( trackIterator - (*trackCollection).begin() );
137 if ( truthIterator == trackTruthCollection->end() ){
138 if (msgLvl(MSG::DEBUG)) msg(MSG::DEBUG) << "No matching truth particle found for track" << endmsg;
139 } else {
140 trackTruth = &((*truthIterator).second);
141 if ( !(trackTruth->particleLink().isValid()) ) {
142 if (msgLvl(MSG::DEBUG)) msg(MSG::DEBUG) << "Link to generated particle information is not there - assuming a lost G4 particle ('fake fake')." << endmsg;
143 // genParticle = m_visibleParticleWithoutTruth; // with pdg_id 0
144 } else {
145#ifdef HEPMC3
146 genParticle = trackTruth->particleLink().scptr();
147#else
148 genParticle = trackTruth->particleLink().cptr();
149#endif
150 if ( genParticle && genParticle->pdg_id() == 0 ) {
151 if (msgLvl(MSG::DEBUG)) msg(MSG::DEBUG) << "Associated Particle ID " << genParticle->pdg_id()
152 << " does not conform to PDG requirements... ignore it!" << endmsg;
153 genParticle = nullptr;
154 }
155 }
156 }
157 if (genParticle) {
158 if (msgLvl(MSG::VERBOSE)) msg(MSG::VERBOSE) <<
159 "Associated Particle ID: " << genParticle->pdg_id() << endmsg;
160 // Perform extrapolation to generate perigee parameters
161 if ( genParticle->production_vertex() )
162 generatedTrackPerigee = m_truthToTrack->makePerigeeParameters( genParticle );
163 }
164
165 // Getting the information if a track has TRT hits or not
166 bool track_has_trthits(false);
167 const Trk::TrackStates *currentTSOSList
168 = (*trackIterator)->trackStateOnSurfaces();
170 = currentTSOSList->begin();
171 for (;itTSOS!=currentTSOSList->end();++itTSOS) {
172 const Trk::MeasurementBase* mb = (*itTSOS)->measurementOnTrack();
173 Identifier id = Trk::IdentifierExtractor::extract(mb);
174 if (id.is_valid() && m_idHelper->is_trt(id)) {
175 track_has_trthits=true;
176 break;
177 }
178 }
179 const TrackParameters* reconstructedPerigee = (*trackIterator)->perigeeParameters();
180
181 if (generatedTrackPerigee != nullptr) {
182
183 if (!reconstructedPerigee) return StatusCode::FAILURE;
184
185 double this_eta = reconstructedPerigee->eta();
186
187 double truth_over_recP = generatedTrackPerigee->parameters()[Trk::qOverP]
188 / reconstructedPerigee->parameters()[Trk::qOverP];
189 if ( (truth_over_recP > 0.9) && (truth_over_recP < 1.1)) {
190 monitorTrackFits(m_nTen,this_eta);
191 monitorTrackFits(m_nTwenty,this_eta);
192 monitorTrackFits(m_nFifty,this_eta);
194 if (track_has_trthits) {
195 monitorTrackFits(m_tTen, this_eta);
196 monitorTrackFits(m_tTwenty, this_eta);
197 monitorTrackFits(m_tFifty, this_eta);
198 monitorTrackFits(m_tHundred, this_eta);
199 }
200
201 } else if ( (truth_over_recP > 0.8) && (truth_over_recP < 1.2) ) {
202 monitorTrackFits(m_nTwenty,this_eta);
203 monitorTrackFits(m_nFifty,this_eta);
205 if (track_has_trthits) {
206 monitorTrackFits(m_tTwenty, this_eta);
207 monitorTrackFits(m_tFifty, this_eta);
208 monitorTrackFits(m_tHundred, this_eta);
209 }
210 } else if ( (truth_over_recP > 0.5) && (truth_over_recP < 1.5) ) {
211 monitorTrackFits(m_nFifty,this_eta);
213 if (track_has_trthits) {
214 monitorTrackFits(m_tFifty, this_eta);
215 monitorTrackFits(m_tHundred, this_eta);
216 }
217 } else {
219 if (track_has_trthits) monitorTrackFits(m_tHundred, this_eta);
220 }
221
222 delete generatedTrackPerigee;
223 } else {
224 double this_eta = reconstructedPerigee->eta();
226 if (track_has_trthits) monitorTrackFits(m_tFakeOrLost, this_eta);
227 }
228
229 } // if it is from a certain track selector/pattern reco info
230
231 } //loop trackCollection end
232
233
234 return StatusCode::SUCCESS;
235}
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
static Double_t sc
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
const HepMcParticleLink & particleLink() const
Definition TrackTruth.h:26
static void extract(std::vector< Identifier > &ids, const std::vector< const MeasurementBase * > &measurements)
static void monitorTrackFits(std::vector< unsigned int > &, const double &)
::StatusCode StatusCode
StatusCode definition for legacy code.
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
DataVector< const Trk::TrackStateOnSurface > TrackStates
@ qOverP
perigee
Definition ParamDefs.h:67
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ 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 Trk::RecMomentumQualityValidation::finalize ( )

standard Athena-Algorithm method

Definition at line 82 of file RecMomentumQualityValidation.cxx.

83{
84 // Code entered here will be executed once at the end of the program run.
85
86 if (msgLvl(MSG::INFO)) printTable();
87
88 return StatusCode::SUCCESS;
89}
void printTable() const
method: make the output table

◆ initialize()

StatusCode Trk::RecMomentumQualityValidation::initialize ( )

standard Athena-Algorithm method

Definition at line 49 of file RecMomentumQualityValidation.cxx.

50{
51
52 StatusCode sc = StatusCode::SUCCESS;
53 if (sc.isFailure()) { /* mute sc*/ }
54
55 // Get the Track Selector Tool
56 if ( !m_trackSelector.empty() ) {
57 sc = m_trackSelector.retrieve();
58 if (sc.isFailure()) {
59 msg(MSG::FATAL) << "Could not retrieve "<< m_trackSelector <<" (to select the tracks which are to be counted) "<< endmsg;
60 msg(MSG::INFO) << "Set the ToolHandle to None if track selection is supposed to be disabled" << endmsg;
61 return sc;
62 }
63 }
64 sc = m_truthToTrack.retrieve();
65 if (sc.isFailure()) {
66 msg(MSG::FATAL) << "Could not retrieve "<< m_truthToTrack << endmsg;
67 return sc;
68 }
69
70 if (detStore()->retrieve(m_idHelper, "AtlasID").isFailure()) {
71 ATH_MSG_ERROR ("Could not get AtlasDetectorID helper" );
72 return StatusCode::FAILURE;
73 }
74
75 ATH_MSG_INFO( "initialize() successful in " << name() );
76 return StatusCode::SUCCESS;
77
78}
#define ATH_MSG_INFO(x)
const ServiceHandle< StoreGateSvc > & detStore() const
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

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

◆ monitorTrackFits()

void Trk::RecMomentumQualityValidation::monitorTrackFits ( std::vector< unsigned int > & Ntracks,
const double & eta )
staticprivate

Definition at line 239 of file RecMomentumQualityValidation.cxx.

240 {
242 if (std::abs(eta) < 0.80 ) ++Ntracks[Trk::RecMomentumQualityValidation::iBarrel];
243 else if (std::abs(eta) < 1.60) ++Ntracks[Trk::RecMomentumQualityValidation::iTransi];
244 else if (std::abs(eta) < 2.40) ++Ntracks[Trk::RecMomentumQualityValidation::iEndcap];
245}
Scalar eta() const
pseudorapidity method

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

◆ printTable()

void Trk::RecMomentumQualityValidation::printTable ( ) const
private

method: make the output table

Definition at line 248 of file RecMomentumQualityValidation.cxx.

248 {
249
250 unsigned int iw = 2;
251 std::string d = " ";
252 for (unsigned int iref=100; iref<m_nHundred[iAll]; iref*=10, ++iw, d+=" ") {}
253 std::cout << "---------------------------------------------------------------------------------" << std::endl;
254 std::cout << " "<< name() << " results "
255 << (m_trackSelector.empty() ? " " : "(with track selection)") << std::endl;
256 std::cout << "---------------------------------------------------------------------------------" << std::endl;
257 std::cout << " q/p truth vicinity -- Any "<<d<<" 50% "<<d<<" 20% "<<d<<" 10%"<<d<<"noTruth " << std::endl;
258 std::cout << "---------------------------------------------------------------------------------" << std::endl;
259 std::cout << " total (Si+TRT) :" << std::setiosflags(std::ios::dec) << std::setw(iw+1)
260 << m_nHundred[iAll] << " ("<< std::setw(iw)<<m_tHundred[iAll]<<") "
261 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
262 << m_nFifty[iAll] << " ("<< std::setw(iw)<<m_tFifty[iAll]<<") "
263 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
264 << m_nTwenty[iAll] << " ("<< std::setw(iw)<<m_tTwenty[iAll]<<") "
265 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
266 << m_nTen[iAll] << " ("<< std::setw(iw)<<m_tTen[iAll]<<") "
267 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
268 << m_nFakeOrLost[iAll] << " ("<< std::setw(iw)<<m_tFakeOrLost[iAll]<<") "
269 << std::endl;
270 std::cout << " barrel (Si+TRT) :" << std::setiosflags(std::ios::dec) << std::setw(iw+1)
271 << m_nHundred[iBarrel] << " ("<< std::setw(iw)<<m_tHundred[iBarrel]<<") "
272 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
273 << m_nFifty[iBarrel] << " ("<< std::setw(iw)<<m_tFifty[iBarrel]<<") "
274 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
275 << m_nTwenty[iBarrel] <<" ("<< std::setw(iw)<<m_tTwenty[iBarrel]<<") "
276 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
277 << m_nTen[iBarrel] << " ("<< std::setw(iw)<<m_tTen[iBarrel]<<") "
278 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
279 << m_nFakeOrLost[iBarrel] << " ("<< std::setw(iw)<<m_tFakeOrLost[iBarrel]<<") "
280 << std::endl;
281 std::cout << " transition :" << std::setiosflags(std::ios::dec) << std::setw(iw+1)
282 << m_nHundred[iTransi] << " ("<< std::setw(iw)<<m_tHundred[iTransi]<<") "
283 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
284 << m_nFifty[iTransi] << " ("<< std::setw(iw)<<m_tFifty[iTransi]<<") "
285 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
286 << m_nTwenty[iTransi] <<" ("<< std::setw(iw)<<m_tTwenty[iTransi]<<") "
287 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
288 << m_nTen[iTransi] << " ("<< std::setw(iw)<<m_tTen[iTransi]<<") "
289 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
290 << m_nFakeOrLost[iTransi] << " ("<< std::setw(iw)<<m_tFakeOrLost[iTransi]<<") "
291 << std::endl;
292 std::cout << " endcap (Si+TRT) :" << std::setiosflags(std::ios::dec) << std::setw(iw+1)
293 << m_nHundred[iEndcap] << " ("<< std::setw(iw)<<m_tHundred[iEndcap]<<") "
294 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
295 << m_nFifty[iEndcap] << " ("<< std::setw(iw)<<m_tFifty[iEndcap]<<") "
296 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
297 << m_nTwenty[iEndcap] <<" ("<< std::setw(iw)<<m_tTwenty[iEndcap]<<") "
298 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
299 << m_nTen[iEndcap] << " ("<< std::setw(iw)<<m_tTen[iEndcap]<<") "
300 << std::setiosflags(std::ios::dec) << std::setw(iw+1)
301 << m_nFakeOrLost[iEndcap] << " ("<< std::setw(iw)<<m_tFakeOrLost[iEndcap]<<") "
302 << std::endl;
303 std::cout << "---------------------------------------------------------------------------------" << std::endl << std::endl;
304}

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

◆ 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}
#define ATH_MSG_WARNING(x)
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_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_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_idHelper

const AtlasDetectorID* Trk::RecMomentumQualityValidation::m_idHelper
private

Definition at line 59 of file RecMomentumQualityValidation.h.

◆ m_inputTrackCollection

std::string Trk::RecMomentumQualityValidation::m_inputTrackCollection
private

properties from JobOptions:

job option: the reco track collection name

Definition at line 55 of file RecMomentumQualityValidation.h.

◆ m_nFakeOrLost

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_nFakeOrLost
private

Definition at line 66 of file RecMomentumQualityValidation.h.

◆ m_nFifty

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_nFifty
private

Definition at line 66 of file RecMomentumQualityValidation.h.

◆ m_nHundred

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_nHundred
private

counters

Definition at line 66 of file RecMomentumQualityValidation.h.

◆ m_nTen

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_nTen
private

Definition at line 66 of file RecMomentumQualityValidation.h.

◆ m_nTwenty

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_nTwenty
private

Definition at line 66 of file RecMomentumQualityValidation.h.

◆ m_tFakeOrLost

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_tFakeOrLost
private

Definition at line 67 of file RecMomentumQualityValidation.h.

◆ m_tFifty

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_tFifty
private

Definition at line 67 of file RecMomentumQualityValidation.h.

◆ m_tHundred

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_tHundred
private

Definition at line 67 of file RecMomentumQualityValidation.h.

◆ m_trackSelector

ToolHandle<Trk::ITrackSelectorTool> Trk::RecMomentumQualityValidation::m_trackSelector
private

Tool handle to Trk::ITrackSelectorTool.

Definition at line 58 of file RecMomentumQualityValidation.h.

◆ m_trackTruthCollection

std::string Trk::RecMomentumQualityValidation::m_trackTruthCollection
private

job option: the truth track collection name

Definition at line 56 of file RecMomentumQualityValidation.h.

◆ m_truthToTrack

ToolHandle<Trk::ITruthToTrack> Trk::RecMomentumQualityValidation::m_truthToTrack
private

Tool handle to Trk::ITruthToTrack tool.

Definition at line 57 of file RecMomentumQualityValidation.h.

◆ m_tTen

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_tTen
private

Definition at line 67 of file RecMomentumQualityValidation.h.

◆ m_tTwenty

std::vector<unsigned int> Trk::RecMomentumQualityValidation::m_tTwenty
private

Definition at line 67 of file RecMomentumQualityValidation.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.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: