ATLAS Offline Software
Loading...
Searching...
No Matches
xAOD::TrackSummary_v1 Class Reference

Track Summary for Acts MultiTrajectory. More...

#include <TrackSummary_v1.h>

Inheritance diagram for xAOD::TrackSummary_v1:
Collaboration diagram for xAOD::TrackSummary_v1:

Public Member Functions

 TrackSummary_v1 ()=default
template<std::size_t measdim = 6>
Eigen::Map< const Eigen::Matrix< double, measdim, 1 > > paramsEigen () const
 access track backend vector of const element
template<std::size_t measdim = 6>
Eigen::Map< Eigen::Matrix< double, measdim, 1 > > paramsEigen ()
 access parameters of non const element
const std::vector< double > & params () const
 access track parameters as plain vector
void setParams (const std::vector< double > &m)
 access set parameters from plain vector
template<std::size_t measdim = 6>
Eigen::Map< const Eigen::Matrix< double, measdim, measdim > > covParamsEigen () const
 access track covariance matrix (flattened, rows layout) of const element
template<std::size_t measdim = 6>
Eigen::Map< Eigen::Matrix< double, measdim, measdim > > covParamsEigen ()
 access track covariance matrix (flattened, rows layout)
const std::vector< double > & covParams () const
 access track covariance as plain vector
void setCovParams (const std::vector< double > &m)
 access set covariance from plain vector
unsigned int nMeasurements () const
 access nMeasurements
void setnMeasurements (unsigned int m)
 access set nMeasurements
const unsigned intnMeasurementsPtr () const
 pointers API needed by MTJ
unsigned intnMeasurementsPtr ()
unsigned int nHoles () const
 access nHoles
void setnHoles (unsigned int m)
 access set nHoles
const unsigned intnHolesPtr () const
 pointers API needed by MTJ
unsigned intnHolesPtr ()
float chi2f () const
 access chi2
void setChi2f (float m)
 access set chi2
const floatchi2fPtr () const
 pointers API needed by MTJ
floatchi2fPtr ()
unsigned int ndf () const
 access ndf
void setNdf (unsigned int m)
 access set ndf
const unsigned intndfPtr () const
 pointers API needed by MTJ
unsigned intndfPtr ()
unsigned int nOutliers () const
 access nOutliers
void setnOutliers (unsigned int m)
 access set nOutliers
const unsigned intnOutliersPtr () const
 pointers API needed by MTJ
unsigned intnOutliersPtr ()
unsigned int nSharedHits () const
 access nSharedHits
void setnSharedHits (unsigned int m)
 access set nSharedHits
const unsigned intnSharedHitsPtr () const
 pointers API needed by MTJ
unsigned intnSharedHitsPtr ()
unsigned int tipIndex () const
 index of the tip of the TrackStates linked list in MultiTrajectory
void setTipIndex (unsigned int)
 Set the tip index.
const unsigned inttipIndexPtr () const
 pointers API needed by MTJ
unsigned inttipIndexPtr ()
unsigned int stemIndex () const
 index of the stem of the TrackStates linked list in MultiTrajectory steemIndex is needed reverse lookup in MultiTrajectory
void setStemIndex (unsigned int)
 Set the stem index.
const unsigned intstemIndexPtr () const
 pointers API needed by MTJ
unsigned intstemIndexPtr ()
unsigned int surfaceIndex () const
 index of the surfaces in the surfaces collection
void setSurfaceIndex (unsigned int)
 Set the index in surface container.
const uint8_tparticleHypothesis () const
 particle hypothesis access
void setParticleHypothesis (const uint8_t &)
void resize (size_t sz=6)
 resize internal arrays to store params (to capacity sz) & convariances (to capacity sz x sz)
size_t size () const
 retrieve the size of the internal vectors for the data summary

Static Private Attributes

static const SG::AuxElement::Accessor< std::vector< double > > s_paramsAcc {"params"}
static const SG::AuxElement::Accessor< std::vector< double > > s_covParamsAcc {"covParams"}

Detailed Description

Track Summary for Acts MultiTrajectory.

Definition at line 20 of file TrackSummary_v1.h.

Constructor & Destructor Documentation

◆ TrackSummary_v1()

xAOD::TrackSummary_v1::TrackSummary_v1 ( )
default

Member Function Documentation

◆ chi2f()

float xAOD::TrackSummary_v1::chi2f ( ) const

access chi2

◆ chi2fPtr() [1/2]

float * xAOD::TrackSummary_v1::chi2fPtr ( )

◆ chi2fPtr() [2/2]

const float * xAOD::TrackSummary_v1::chi2fPtr ( ) const

pointers API needed by MTJ

◆ covParams()

const std::vector< double > & xAOD::TrackSummary_v1::covParams ( ) const

access track covariance as plain vector

◆ covParamsEigen() [1/2]

template<std::size_t measdim = 6>
Eigen::Map< Eigen::Matrix< double, measdim, measdim > > xAOD::TrackSummary_v1::covParamsEigen ( )
inline

access track covariance matrix (flattened, rows layout)

Definition at line 66 of file TrackSummary_v1.h.

67 {
68 return Eigen::Map<Eigen::Matrix<double, measdim, measdim>>{s_covParamsAcc(*this).data()};
69 }
static const SG::AuxElement::Accessor< std::vector< double > > s_covParamsAcc

◆ covParamsEigen() [2/2]

template<std::size_t measdim = 6>
Eigen::Map< const Eigen::Matrix< double, measdim, measdim > > xAOD::TrackSummary_v1::covParamsEigen ( ) const
inline

access track covariance matrix (flattened, rows layout) of const element

Definition at line 57 of file TrackSummary_v1.h.

58 {
59 return Eigen::Map<const Eigen::Matrix<double, measdim, measdim>>{s_covParamsAcc(*this).data()};
60 }

◆ ndf()

unsigned int xAOD::TrackSummary_v1::ndf ( ) const

access ndf

◆ ndfPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::ndfPtr ( )

◆ ndfPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::ndfPtr ( ) const

pointers API needed by MTJ

◆ nHoles()

unsigned int xAOD::TrackSummary_v1::nHoles ( ) const

access nHoles

◆ nHolesPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::nHolesPtr ( )

◆ nHolesPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::nHolesPtr ( ) const

pointers API needed by MTJ

◆ nMeasurements()

unsigned int xAOD::TrackSummary_v1::nMeasurements ( ) const

access nMeasurements

◆ nMeasurementsPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::nMeasurementsPtr ( )

◆ nMeasurementsPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::nMeasurementsPtr ( ) const

pointers API needed by MTJ

◆ nOutliers()

unsigned int xAOD::TrackSummary_v1::nOutliers ( ) const

access nOutliers

◆ nOutliersPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::nOutliersPtr ( )

◆ nOutliersPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::nOutliersPtr ( ) const

pointers API needed by MTJ

◆ nSharedHits()

unsigned int xAOD::TrackSummary_v1::nSharedHits ( ) const

access nSharedHits

◆ nSharedHitsPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::nSharedHitsPtr ( )

◆ nSharedHitsPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::nSharedHitsPtr ( ) const

pointers API needed by MTJ

◆ params()

const std::vector< double > & xAOD::TrackSummary_v1::params ( ) const

access track parameters as plain vector

◆ paramsEigen() [1/2]

template<std::size_t measdim = 6>
Eigen::Map< Eigen::Matrix< double, measdim, 1 > > xAOD::TrackSummary_v1::paramsEigen ( )
inline

access parameters of non const element

Definition at line 39 of file TrackSummary_v1.h.

40 {
41 return Eigen::Map<Eigen::Matrix<double, measdim, 1>>{s_paramsAcc(*this).data()};
42 }
static const SG::AuxElement::Accessor< std::vector< double > > s_paramsAcc

◆ paramsEigen() [2/2]

template<std::size_t measdim = 6>
Eigen::Map< const Eigen::Matrix< double, measdim, 1 > > xAOD::TrackSummary_v1::paramsEigen ( ) const
inline

access track backend vector of const element

Definition at line 31 of file TrackSummary_v1.h.

32 {
33 return Eigen::Map<const Eigen::Matrix<double, measdim, 1>>{s_paramsAcc(*this).data()};
34 }

◆ particleHypothesis()

const uint8_t & xAOD::TrackSummary_v1::particleHypothesis ( ) const

particle hypothesis access

◆ resize()

void xAOD::TrackSummary_v1::resize ( size_t sz = 6)

resize internal arrays to store params (to capacity sz) & convariances (to capacity sz x sz)

Definition at line 45 of file TrackSummary_v1.cxx.

45 {
46 s_paramsAcc(*this).resize(sz);
47 s_covParamsAcc(*this).resize(sz * sz);
48}
static Double_t sz

◆ setChi2f()

void xAOD::TrackSummary_v1::setChi2f ( float m)

access set chi2

◆ setCovParams()

void xAOD::TrackSummary_v1::setCovParams ( const std::vector< double > & m)

access set covariance from plain vector

◆ setNdf()

void xAOD::TrackSummary_v1::setNdf ( unsigned int m)

access set ndf

◆ setnHoles()

void xAOD::TrackSummary_v1::setnHoles ( unsigned int m)

access set nHoles

◆ setnMeasurements()

void xAOD::TrackSummary_v1::setnMeasurements ( unsigned int m)

access set nMeasurements

◆ setnOutliers()

void xAOD::TrackSummary_v1::setnOutliers ( unsigned int m)

access set nOutliers

◆ setnSharedHits()

void xAOD::TrackSummary_v1::setnSharedHits ( unsigned int m)

access set nSharedHits

◆ setParams()

void xAOD::TrackSummary_v1::setParams ( const std::vector< double > & m)

access set parameters from plain vector

◆ setParticleHypothesis()

void xAOD::TrackSummary_v1::setParticleHypothesis ( const uint8_t & )

◆ setStemIndex()

void xAOD::TrackSummary_v1::setStemIndex ( unsigned int )

Set the stem index.

See also
stemIndex() for explanation

◆ setSurfaceIndex()

void xAOD::TrackSummary_v1::setSurfaceIndex ( unsigned int )

Set the index in surface container.

◆ setTipIndex()

void xAOD::TrackSummary_v1::setTipIndex ( unsigned int )

Set the tip index.

See also
tipIndex() for explanation

◆ size()

size_t xAOD::TrackSummary_v1::size ( ) const

retrieve the size of the internal vectors for the data summary

Definition at line 50 of file TrackSummary_v1.cxx.

50 {
51 return s_paramsAcc(*this).size();
52}

◆ stemIndex()

unsigned int xAOD::TrackSummary_v1::stemIndex ( ) const

index of the stem of the TrackStates linked list in MultiTrajectory steemIndex is needed reverse lookup in MultiTrajectory

◆ stemIndexPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::stemIndexPtr ( )

◆ stemIndexPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::stemIndexPtr ( ) const

pointers API needed by MTJ

◆ surfaceIndex()

unsigned int xAOD::TrackSummary_v1::surfaceIndex ( ) const

index of the surfaces in the surfaces collection

◆ tipIndex()

unsigned int xAOD::TrackSummary_v1::tipIndex ( ) const

index of the tip of the TrackStates linked list in MultiTrajectory

◆ tipIndexPtr() [1/2]

unsigned int * xAOD::TrackSummary_v1::tipIndexPtr ( )

◆ tipIndexPtr() [2/2]

const unsigned int * xAOD::TrackSummary_v1::tipIndexPtr ( ) const

pointers API needed by MTJ

Member Data Documentation

◆ s_covParamsAcc

const SG::AuxElement::Accessor< std::vector< double > > xAOD::TrackSummary_v1::s_covParamsAcc {"covParams"}
staticprivate

Definition at line 43 of file TrackSummary_v1.h.

◆ s_paramsAcc

setCovParams const SG::AuxElement::Accessor< std::vector< double > > xAOD::TrackSummary_v1::s_paramsAcc {"params"}
staticprivate

Definition at line 41 of file TrackSummary_v1.h.

41{s_paramsAcc(*this).data()};

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