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

This validation tool writes information about Trk::CompetingRIOsOnTrack into an ntuple. More...

#include <DAF_ValidationNtupleHelper.h>

Inheritance diagram for Trk::DAF_ValidationNtupleHelper:
Collaboration diagram for Trk::DAF_ValidationNtupleHelper:

Public Member Functions

 DAF_ValidationNtupleHelper (const std::string &, const std::string &, const IInterface *)
 ~DAF_ValidationNtupleHelper ()
StatusCode initialize ()
 initialize
StatusCode finalize ()
 finalize
virtual StatusCode fillMeasurementData (const Trk::MeasurementBase *, const Trk::TrackParameters *, const int &detectorType, const bool &isOutlier)
 fill Trk::CompetingRIOsOnTrack data
virtual StatusCode fillHoleData (const Trk::TrackStateOnSurface &, const int &)
 fill special data about holes on track (here: do nothing)
virtual StatusCode addNtupleItems (TTree *tree, const int &detectorType)
 add items to the ntuple and configure the helper tool: should be called once (per detector type) by the steering tool (Trk::ITrackValidationNtupleTool)
virtual StatusCode resetVariables (const int &detectorType)
 reset ntuple variables
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 ()
 Interface ID, declared here, and defined below.

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

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

Private Attributes

bool m_ignoreMissTrkCov
bool m_writeHitPositions
 jobOption: shall the positions of the contained ROTs be written?
ToolHandle< Trk::IResidualPullCalculatorm_residualPullCalculator
 The residual and pull calculator tool.
int * m_isUnbiased
std::vector< int > * m_nContainedROTs
std::vector< int > * m_indexOfMaxAssgnProb
std::vector< float > * m_maxAssgnProb
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

Static Private Attributes

static const unsigned int s_maxContainedROTs = 8

Detailed Description

This validation tool writes information about Trk::CompetingRIOsOnTrack into an ntuple.

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

◆ DAF_ValidationNtupleHelper()

Trk::DAF_ValidationNtupleHelper::DAF_ValidationNtupleHelper ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 35 of file DAF_ValidationNtupleHelper.cxx.

39 :
40AthAlgTool(t,n,p),
41 m_residualPullCalculator("Trk::ResidualPullCalculator/ResidualPullCalculator"),
42 m_isUnbiased(nullptr),
43 m_nContainedROTs(nullptr),
44 m_indexOfMaxAssgnProb(nullptr),
45 m_maxAssgnProb(nullptr)
46
47
48 {
49 declareInterface<IValidationNtupleHelperTool>(this);
50 declareProperty("IgnoreMissingTrackCovarianceForPulls", m_ignoreMissTrkCov = false, "Do not warn, if track states do not have covariance matries when calculating pulls");
51 declareProperty("WriteMeasurementPositionOfROTs", m_writeHitPositions = true, "Write measurement positions?");
52 declareProperty("ResidualPullCalculatorTool", m_residualPullCalculator, "Tool to calculate residuals and pulls");
53}
AthAlgTool()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ToolHandle< Trk::IResidualPullCalculator > m_residualPullCalculator
The residual and pull calculator tool.
bool m_writeHitPositions
jobOption: shall the positions of the contained ROTs be written?

◆ ~DAF_ValidationNtupleHelper()

Trk::DAF_ValidationNtupleHelper::~DAF_ValidationNtupleHelper ( )
default

Member Function Documentation

◆ addNtupleItems()

StatusCode Trk::DAF_ValidationNtupleHelper::addNtupleItems ( TTree * tree,
const int & detectorType )
virtual

add items to the ntuple and configure the helper tool: should be called once (per detector type) by the steering tool (Trk::ITrackValidationNtupleTool)

addNtupleItems

Implements Trk::IValidationNtupleHelperTool.

Definition at line 89 of file DAF_ValidationNtupleHelper.cxx.

92 {
93
94 // retrieve item (is the used track state unbiased?) from ntuple
95 TBranch* trackStatesUnbiasedBranch = tree->GetBranch("TrackStatesUnbiased");
96 if (!trackStatesUnbiasedBranch) {
97 ATH_MSG_ERROR ("Unable to get Branch TrackStatesUnbiased in ntuple" );
98 return StatusCode::FAILURE;
99 }
100 void* variableAdr = static_cast<void*>(trackStatesUnbiasedBranch->GetAddress());
101 if (!variableAdr) {
102 ATH_MSG_ERROR ( "Unable to get variable address of Branch TrackStatesUnbiased" );
103 return StatusCode::FAILURE;
104 }
105 ATH_MSG_VERBOSE ("Address of variable: " << variableAdr);
106 m_isUnbiased = static_cast<int*>(variableAdr);
107
108 // add items
109 m_nContainedROTs = new std::vector<int>();
110 m_indexOfMaxAssgnProb = new std::vector<int>();
111 m_maxAssgnProb = new std::vector<float>();
112 tree->Branch("CompROTnContainedROTs", &m_nContainedROTs); // number of contained ROTs in the CompetingRIOsOnTrack (1-dim array: CompROTnContainedROTs[nHits])
113 tree->Branch("CompROTindexOfMaxAssgnProb", &m_indexOfMaxAssgnProb); // index of the ROT with the maximum assignment probability
114 tree->Branch("CompROTmaxAssgnProb", &m_maxAssgnProb); // maximum assignment probability
115
116 ATH_MSG_VERBOSE ("added items to ntuple");
117 return StatusCode::SUCCESS;
118}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
TChain * tree

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

◆ 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

◆ fillHoleData()

StatusCode Trk::DAF_ValidationNtupleHelper::fillHoleData ( const Trk::TrackStateOnSurface & ,
const int &  )
virtual

fill special data about holes on track (here: do nothing)

Implements Trk::IValidationNtupleHelperTool.

Definition at line 167 of file DAF_ValidationNtupleHelper.cxx.

169 {
170 // we do nothing with holes
171 return StatusCode::SUCCESS;
172}

◆ fillMeasurementData()

StatusCode Trk::DAF_ValidationNtupleHelper::fillMeasurementData ( const Trk::MeasurementBase * measurement,
const Trk::TrackParameters * trkParameters,
const int & detectorType,
const bool & isOutlier )
virtual

fill Trk::CompetingRIOsOnTrack data

fill CompetingRIOsOnTrack data

Implements Trk::IValidationNtupleHelperTool.

Definition at line 123 of file DAF_ValidationNtupleHelper.cxx.

129 {
130
131 if (isOutlier) {
132 // we do nothing for outliers, because competing ROTs are no outliers by definition
133 return StatusCode::SUCCESS;
134 }
135 if (!trkParameters) {
136 ATH_MSG_ERROR ( "no TrackParameters given" );
137 return StatusCode::FAILURE;
138 }
139
140 const Trk::CompetingRIOsOnTrack* comprot = dynamic_cast<const Trk::CompetingRIOsOnTrack*>(measurement);
141 if (comprot) {
142 ATH_MSG_VERBOSE ( "measurement is CompetingRIOsOnTrack, fill ntuple variables" );
143 //--------------------------
144 // status info about compROT
145 m_nContainedROTs->push_back(comprot->numberOfContainedROTs());
146 unsigned int indexOfMaxAssgnProb = comprot->indexOfMaxAssignProb();
147 m_indexOfMaxAssgnProb->push_back(indexOfMaxAssgnProb);
148 m_maxAssgnProb->push_back(comprot->assignmentProbability(indexOfMaxAssgnProb));
149
150 } else {
151 m_nContainedROTs->push_back(0);
152 m_indexOfMaxAssgnProb->push_back(0);
153 m_maxAssgnProb->push_back(0.);
154 ATH_MSG_VERBOSE ( "measurement is not a CompetingRIOsOnTrack, ignore it");
155 }
156 return StatusCode::SUCCESS;
157}
AssignmentProb assignmentProbability(unsigned int indx) const
returns the AssignmentProbability depending on the integer.
virtual unsigned int numberOfContainedROTs() const =0
Number of RIO_OnTracks to be contained by this CompetingRIOsOnTrack.
unsigned int indexOfMaxAssignProb() const
Index of the ROT with the highest assignment probability.

◆ finalize()

StatusCode Trk::DAF_ValidationNtupleHelper::finalize ( )

finalize

Definition at line 80 of file DAF_ValidationNtupleHelper.cxx.

80 {
81 ATH_MSG_INFO ( "finalize() successful in " << name() );
82 return StatusCode::SUCCESS;
83}
#define ATH_MSG_INFO(x)

◆ initialize()

StatusCode Trk::DAF_ValidationNtupleHelper::initialize ( )

initialize

Definition at line 63 of file DAF_ValidationNtupleHelper.cxx.

63 {
64
66
67 // get the ResidualPullCalculator
68 sc = m_residualPullCalculator.retrieve();
69 if (sc.isFailure()) {
70 ATH_MSG_FATAL ("Could not retrieve "<< m_residualPullCalculator <<" (to calculate residuals and pulls) " );
71 return sc;
72 }
73
74 return StatusCode::SUCCESS;
75}
#define ATH_MSG_FATAL(x)
static Double_t sc
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ 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 & Trk::IValidationNtupleHelperTool::interfaceID ( )
inlinestaticinherited

Interface ID, declared here, and defined below.

Definition at line 76 of file IValidationNtupleHelperTool.h.

76 {
78}
static const InterfaceID IID_IValidationNtupleHelperTool("IValidationNtupleHelperTool", 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 }

◆ resetVariables()

StatusCode Trk::DAF_ValidationNtupleHelper::resetVariables ( const int & detectorType)
virtual

reset ntuple variables

Implements Trk::IValidationNtupleHelperTool.

Definition at line 159 of file DAF_ValidationNtupleHelper.cxx.

160 {
161 m_nContainedROTs->clear();
162 m_indexOfMaxAssgnProb->clear();
163 m_maxAssgnProb->clear();
164 return StatusCode::SUCCESS;
165}

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

bool Trk::DAF_ValidationNtupleHelper::m_ignoreMissTrkCov
private

Definition at line 72 of file DAF_ValidationNtupleHelper.h.

◆ m_indexOfMaxAssgnProb

std::vector<int>* Trk::DAF_ValidationNtupleHelper::m_indexOfMaxAssgnProb
private

Definition at line 80 of file DAF_ValidationNtupleHelper.h.

◆ m_isUnbiased

int* Trk::DAF_ValidationNtupleHelper::m_isUnbiased
private

Definition at line 77 of file DAF_ValidationNtupleHelper.h.

◆ m_maxAssgnProb

std::vector<float>* Trk::DAF_ValidationNtupleHelper::m_maxAssgnProb
private

Definition at line 81 of file DAF_ValidationNtupleHelper.h.

◆ m_nContainedROTs

std::vector<int>* Trk::DAF_ValidationNtupleHelper::m_nContainedROTs
private

Definition at line 79 of file DAF_ValidationNtupleHelper.h.

◆ m_residualPullCalculator

ToolHandle<Trk::IResidualPullCalculator> Trk::DAF_ValidationNtupleHelper::m_residualPullCalculator
private

The residual and pull calculator tool.

Definition at line 74 of file DAF_ValidationNtupleHelper.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_writeHitPositions

bool Trk::DAF_ValidationNtupleHelper::m_writeHitPositions
private

jobOption: shall the positions of the contained ROTs be written?

Definition at line 73 of file DAF_ValidationNtupleHelper.h.

◆ s_maxContainedROTs

const unsigned int Trk::DAF_ValidationNtupleHelper::s_maxContainedROTs = 8
staticprivate

Definition at line 71 of file DAF_ValidationNtupleHelper.h.


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