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

Class describing a TrackParticle. More...

#include <TrackParticle_v1.h>

Inheritance diagram for xAOD::TrackParticle_v1:
Collaboration diagram for xAOD::TrackParticle_v1:

Public Member Functions

 TrackParticle_v1 ()
 Default constructor.
 ~TrackParticle_v1 ()
 Destructor.
 TrackParticle_v1 (const TrackParticle_v1 &o)
 Copy ctor. This involves copying the entire Auxilary store, and is a slow operation which should be used sparingly.
TrackParticle_v1operator= (const TrackParticle_v1 &tp)
 Assignment operator. This can involve creating and copying an Auxilary store, and so should be used sparingly.
bool summaryValue (uint8_t &value, const SummaryType &information) const
 Accessor for TrackSummary values.
bool summaryValue (float &value, const SummaryType &information) const
 Accessor for TrackSummary values.
void setSummaryValue (uint8_t &value, const SummaryType &information)
 Set method for TrackSummary values.
void setSummaryValue (float &value, const SummaryType &information)
 Set method for TrackSummary values.
void resetCache ()
 Reset the internal cache of the object.

Private Types

using covMatrixIndex = TrackingDetails::covMatrixIndex
typedef std::vector< std::pair< covMatrixIndex, covMatrixIndex > > covMatrixIndexPairVec

Static Private Member Functions

static const covMatrixIndexPairVeccovMatrixComprIndexPairs ()

Private Attributes

CxxUtils::CachedValue< Trk::Perigeem_perigeeParameters
 Cached MeasuredPerigee, built from this object.

Static Private Attributes

static const std::size_t COVMATRIX_OFFDIAG_VEC_COMPR_SIZE = TrackingDetails::COVMATRIX_OFFDIAG_VEC_COMPR_SIZE

IParticle functions

typedef IParticle::FourMom_t FourMom_t
 Definition of the 4-momentum type.
typedef ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D< double > > GenVecFourMom_t
 Base 4 Momentum type for TrackParticle.
virtual double pt () const override final
 The transverse momentum ( \(p_T\)) of the particle.
virtual double eta () const override final
 The pseudorapidity ( \(\eta\)) of the particle.
virtual double phi () const override final
 The azimuthal angle ( \(\phi\)) of the particle (has range \(-\pi\) to \(+\pi\).).
virtual double m () const override final
 The invariant mass of the particle..
virtual double e () const override final
 The total energy of the particle.
virtual double rapidity () const override final
 The true rapidity (y) of the particle.
virtual FourMom_t p4 () const override final
 The full 4-momentum of the particle.
GenVecFourMom_t genvecP4 () const
 The full 4-momentum of the particle : GenVector form.
virtual Type::ObjectType type () const override final
 The type of the object as a simple enumeration.

Defining parameters functions

The 'defining parameters' are key to the concept of a TrackParticle, and give the values for the IParticle interface ( pt(), phi(), eta() etc.).

They use the Trk::Perigee coordinate system, and are defined as: \(( d_0, z_0, \phi, \theta, q/p )\). The parameters are expressed with respect to an origin (returned by vx(), vy() and vy() ), currently intended to be the 'beamspot'. This origin is expected to be the same for all track particles in a collection (and this may be be enforced). The \(\phi\) parameter is returned by either the phi() or the phi0() methods, the difference just being whether it is returned as a float or a double (it is stored as a float)

float charge () const
 Returns the charge.
float d0 () const
 Returns the \(d_0\) parameter.
float z0 () const
 Returns the \(z_0\) parameter.
float phi0 () const
 Returns the \(\phi\) parameter, which has range \(-\pi\) to \(+\pi\).
float theta () const
 Returns the \(\theta\) parameter, which has range 0 to \(\pi\).
float qOverP () const
 Returns the \(q/p\) parameter.
float time () const
 Returns the time.
float timeResolution () const
 Returns the time resolution.
uint8_t hasValidTime () const
 Returns whether or not the track has a valid time.
DefiningParameters_t definingParameters () const
 Returns a SVector of the Perigee track parameters.
const ParametersCovMatrix_t definingParametersCovMatrix () const
 Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.
ParametersCovMatrixFilled_t definingParametersCovMatrixFilled () const
 Returns a 5x5 matrix describing which elements of the covariance matrix are known.
const std::vector< float > & definingParametersCovMatrixDiagVec () const
 Returns the diagonal elements of the defining parameters covariance matrix.
const std::vector< float > & definingParametersCovMatrixOffDiagVec () const
 Returns the correlation coefficient associated with the off-diagonal elements of the covariance matrix = cov(X,Y)/sqrt(cov(X,X)*cov(Y,Y)).
std::vector< floatdefiningParametersCovMatrixVec () const
 Returns the length 6 vector containing the elements of defining parameters covariance matrix.
bool definingParametersCovMatrixOffDiagCompr () const
void setDefiningParameters (float d0, float z0, float phi0, float theta, float qOverP)
 Set the defining parameters.
void setDefiningParameters (float d0, float z0, float phi0, float theta, float qOverP, float time)
void setTime (float time)
void setTimeResolution (float timeResolution)
void setHasValidTime (uint8_t hasValidTime)
void setDefiningParametersCovMatrix (const ParametersCovMatrix_t &cov)
 Set the defining parameters covariance matrix.
void setDefiningParametersCovMatrixDiagVec (const std::vector< float > &vec)
 Set the defining parameters covariance matrix using a length 15 vector.
void setDefiningParametersCovMatrixOffDiagVec (const std::vector< float > &vec)
 Set the off-diagonal elements of the defining parameters covariance matrix.
void setDefiningParametersCovMatrixVec (const std::vector< float > &cov)
void compressDefiningParametersCovMatrixOffDiag ()
 Delete some off-diagonal elements for compression.
float vx () const
 The x origin for the parameters.
float vy () const
 The y origin for the parameters.
float vz () const
 The z origin for the parameters.
void setParametersOrigin (float x, float y, float z)
 Set the origin for the parameters.
const Trk::PerigeeperigeeParameters () const
 Returns the Trk::MeasuredPerigee track parameters.

Curvilinear functions

The set of functions which return other track parameters.

The remaining track parameters (i.e. not the 'defining parameters') use the 'curvilinear' coordinate system, and are represented by the parameters \((x,y,z,p_x,p_y,p_z)\). The parameters can have an associated local 5x5 error/covariance matrix. They are expressed at various points through the detector, which can be determined by the parameterPosition() method.

// Example code to use parameters
unsigned int index=0;
if (myTP.indexOfParameterAtPosition(index, xAOD::FirstMeasurement)){
CurvilinearParameters_t parameters = myTP.trackParameters(index);
}
Definition index.py:1
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.
size_t numberOfParameters () const
 Returns the number of additional parameters stored in the TrackParticle.
const CurvilinearParameters_t trackParameters (unsigned int index) const
 Returns the track parameter vector at 'index'.
float parameterX (unsigned int index) const
 Returns the parameter x position, for 'index'.
float parameterY (unsigned int index) const
 Returns the parameter y position, for 'index'.
float parameterZ (unsigned int index) const
 Returns the parameter z position, for 'index'.
float parameterPX (unsigned int index) const
 Returns the parameter x momentum component, for 'index'.
float parameterPY (unsigned int index) const
 Returns the parameter y momentum component, for 'index'.
float parameterPZ (unsigned int index) const
 Returns the parameter z momentum component, for 'index'.
void setTrackParameters (std::vector< std::vector< float > > &parameters)
 Set the parameters via the passed vector of vectors.
ParametersCovMatrix_t trackParameterCovarianceMatrix (unsigned int index) const
 Returns the TrackParticleCovMatrix_t (covariance matrix) at 'index', which corresponds to the parameters at the same index.
void setTrackParameterCovarianceMatrix (unsigned int index, std::vector< float > &cov)
 Set the cov matrix of the parameter at 'index', using a vector of floats.
xAOD::ParameterPosition parameterPosition (unsigned int index) const
 Return the ParameterPosition of the parameters at 'index'.
bool indexOfParameterAtPosition (unsigned int &index, ParameterPosition position) const
 Function to determine if this TrackParticle contains track parameters at a certain position, and if so, what the 'index' is.
void setParameterPosition (unsigned int index, ParameterPosition pos)
 Set the 'position' (i.e. where it is in ATLAS) of the parameter at 'index', using the ParameterPosition enum.
const Trk::CurvilinearParameters curvilinearParameters (unsigned int index) const
 Returns a curvilinear representation of the parameters at 'index'.
float radiusOfFirstHit () const
 Returns the radius of the first hit.
void setRadiusOfFirstHit (float radius)
 Set the radius of the first hit.
uint64_t identifierOfFirstHit () const
 Returns the offline identifier of the first hit.
void setIdentifierOfFirstHit (uint64_t id)
 Set the offline identifier of the first hit.
float beamlineTiltX () const
void setBeamlineTiltX (float tiltX)
float beamlineTiltY () const
void setBeamlineTiltY (float tiltY)
uint32_t hitPattern () const
void setHitPattern (uint32_t hitpattern)
uint8_t numberOfUsedHitsdEdx () const
void setNumberOfUsedHitsdEdx (uint8_t numhits)
uint8_t numberOfIBLOverflowsdEdx () const
void setNumberOfIBLOverflowsdEdx (uint8_t numoverflows)

Fit quality functions

Returns some information about quality of the track fit.

float chiSquared () const
 Returns the \( \chi^2 \) of the overall track fit.
float numberDoF () const
 Returns the number of degrees of freedom of the overall track or vertex fit as float.
void setFitQuality (float chiSquared, float numberDoF)
 Set the 'Fit Quality' information.

TrackInfo functions

Contains information about the 'fitter' of this Trk::Track / TrackParticle.

Additionally there is some information about how the e.g. fit was configured. Also the information on the properties of the track fit is stored.

void setTrackProperties (const TrackProperties properties)
 Methods setting the TrackProperties.
void setPatternRecognitionInfo (const std::bitset< xAOD::NumberOfTrackRecoInfo > &patternReco)
 Method setting the pattern recognition algorithm, using a bitset.
void setPatternRecognitionInfo (uint64_t patternReco)
 Method setting the pattern recognition algorithm, using a 64-bit int (which is faster than using a bitset).
void setTrackFitter (const TrackFitter fitter)
 Method for setting the fitter, using the TrackFitter enum.
void setParticleHypothesis (const ParticleHypothesis hypo)
 Method for setting the particle type, using the ParticleHypothesis enum.
TrackProperties trackProperties () const
 Access methods for track properties, which returns 'true' if a logical AND of the parameter 'proprty' and the stored properties returns true.
std::bitset< NumberOfTrackRecoInfopatternRecoInfo () const
 Access method for pattern recognition algorithm.
ParticleHypothesis particleHypothesis () const
 Returns the particle hypothesis used for Track fitting.
TrackFitter trackFitter () const
 Returns the fitter.

Links

const ElementLink< TrackCollection > & trackLink () const
 Returns a link (which can be invalid) to the Trk::Track which was used to make this TrackParticle.
void setTrackLink (const ElementLink< TrackCollection > &track)
 Set the link to the original track.
const Trk::Tracktrack () const
 Returns a pointer (which can be NULL) to the Trk::Track which was used to make this TrackParticle.

Functions for getting and setting user properties

template<class T>
XAOD_AUXDATA_DEPRECATED T & auxdata (const std::string &name, const std::string &clsname="")
 Fetch an aux data variable, as a non-const reference.
template<class T>
XAOD_AUXDATA_DEPRECATED const T & auxdata (const std::string &name, const std::string &clsname="") const
 Fetch an aux data variable, as a const reference.
template<class T>
XAOD_AUXDATA_DEPRECATED bool isAvailable (const std::string &name, const std::string &clsname="") const
 Check if a user property is available for reading or not.
template<class T>
XAOD_AUXDATA_DEPRECATED bool isAvailableWritable (const std::string &name, const std::string &clsname="") const
 Check if a user property is available for writing or not.

Detailed Description

Member Typedef Documentation

◆ covMatrixIndex

◆ covMatrixIndexPairVec

Definition at line 351 of file TrackParticle_v1.h.

◆ FourMom_t

Definition of the 4-momentum type.

Definition at line 73 of file TrackParticle_v1.h.

◆ GenVecFourMom_t

typedef ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > xAOD::TrackParticle_v1::GenVecFourMom_t

Base 4 Momentum type for TrackParticle.

Definition at line 79 of file TrackParticle_v1.h.

Constructor & Destructor Documentation

◆ TrackParticle_v1() [1/2]

xAOD::TrackParticle_v1::TrackParticle_v1 ( )

Default constructor.

Definition at line 47 of file TrackParticle_v1.cxx.

48 : IParticle() {
49 // perigeeParameters cache initialized to be empty (default constructor)
50 }
IParticle()=default

◆ ~TrackParticle_v1()

xAOD::TrackParticle_v1::~TrackParticle_v1 ( )

Destructor.

Definition at line 73 of file TrackParticle_v1.cxx.

73{}

◆ TrackParticle_v1() [2/2]

xAOD::TrackParticle_v1::TrackParticle_v1 ( const TrackParticle_v1 & o)

Copy ctor. This involves copying the entire Auxilary store, and is a slow operation which should be used sparingly.

Definition at line 52 of file TrackParticle_v1.cxx.

53 : IParticle( tp ) {
54 makePrivateStore( tp );
55 // perigeeParameters cache initialized to be empty (default constructor)
56 // assume that this copy will create new cache as needed
57 }
void makePrivateStore()
Create a new (empty) private store for this object.

Member Function Documentation

◆ auxdata() [1/2]

template<class T>
XAOD_AUXDATA_DEPRECATED T & xAOD::IParticle::auxdata ( const std::string & name,
const std::string & clsname = "" )
inlineinherited

Fetch an aux data variable, as a non-const reference.

This function provides an easy way for users to decorate objects with auxiliary data.

Take note that this function is slow. Should not be used inside time-critical code.

Parameters
nameName of the aux variable
clsnameThe name of the associated class. May be blank
Returns
A modifiable reference to the decoration

Definition at line 98 of file Event/xAOD/xAODBase/xAODBase/IParticle.h.

99 {
100
101 return SG::Accessor< T >(name, clsname)(*this);
102 }
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:573

◆ auxdata() [2/2]

template<class T>
XAOD_AUXDATA_DEPRECATED const T & xAOD::IParticle::auxdata ( const std::string & name,
const std::string & clsname = "" ) const
inlineinherited

Fetch an aux data variable, as a const reference.

This function provides an easy way for users to retrieve auxiliary decorations from an object.

Take note that this function is slow. Should not be used inside time-critical code.

Parameters
nameName of the aux variable
clsnameThe name of the associated class. May be blank
Returns
A constant reference to the decoration

Definition at line 118 of file Event/xAOD/xAODBase/xAODBase/IParticle.h.

119 {
120
121 return SG::ConstAccessor< T >( name, clsname )( *this );
122 }
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:570

◆ beamlineTiltX()

float xAOD::TrackParticle_v1::beamlineTiltX ( ) const

◆ beamlineTiltY()

float xAOD::TrackParticle_v1::beamlineTiltY ( ) const

◆ charge()

float xAOD::TrackParticle_v1::charge ( ) const

Returns the charge.

Definition at line 143 of file TrackParticle_v1.cxx.

143 {
145 }
float qOverP() const
Returns the parameter.
float charge(float qOverP)

◆ chiSquared()

float xAOD::TrackParticle_v1::chiSquared ( ) const

Returns the \( \chi^2 \) of the overall track fit.

◆ compressDefiningParametersCovMatrixOffDiag()

void xAOD::TrackParticle_v1::compressDefiningParametersCovMatrixOffDiag ( )

Delete some off-diagonal elements for compression.

Definition at line 388 of file TrackParticle_v1.cxx.

388 {
389
390 ParametersCovMatrix_t cov = definingParametersCovMatrix();
391 std::vector< float > offDiagVecCompr;
392 offDiagVecCompr.resize(COVMATRIX_OFFDIAG_VEC_COMPR_SIZE);
393
394 const covMatrixIndexPairVec& vecPairIndex = covMatrixComprIndexPairs();
395
396 for(unsigned int k=0; k<COVMATRIX_OFFDIAG_VEC_COMPR_SIZE; ++k){
397 std::pair<covMatrixIndex,covMatrixIndex> pairIndex = vecPairIndex[k];
398 covMatrixIndex i = pairIndex.first;
399 covMatrixIndex j = pairIndex.second;
400 float offDiagElement = cov(i,i)>0 && cov(j,j)>0 ? cov(i,j)/sqrt(cov(i,i)*cov(j,j)) : 0;
401 offDiagVecCompr[k] = offDiagElement;
402 }
403
404 accCovMatrixOffDiag( *this ) = offDiagVecCompr;
405 return;
406
407 }
TrackingDetails::covMatrixIndex covMatrixIndex
static const covMatrixIndexPairVec & covMatrixComprIndexPairs()
const ParametersCovMatrix_t definingParametersCovMatrix() const
Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.
static const std::size_t COVMATRIX_OFFDIAG_VEC_COMPR_SIZE
std::vector< std::pair< covMatrixIndex, covMatrixIndex > > covMatrixIndexPairVec
float j(const xAOD::IParticle &, const xAOD::TrackMeasurementValidation &hit, const Eigen::Matrix3d &jab_inv)
static const SG::AuxElement::Accessor< std::vector< float > > accCovMatrixOffDiag("definingParametersCovMatrixOffDiag")

◆ covMatrixComprIndexPairs()

const TrackParticle_v1::covMatrixIndexPairVec & xAOD::TrackParticle_v1::covMatrixComprIndexPairs ( )
staticprivate

Definition at line 716 of file TrackParticle_v1.cxx.

716 {
719 return result;
720 }
constexpr std::array< std::pair< covMatrixIndex, covMatrixIndex >, COVMATRIX_OFFDIAG_VEC_COMPR_SIZE > covMatrixComprIndexPairs

◆ curvilinearParameters()

const Trk::CurvilinearParameters xAOD::TrackParticle_v1::curvilinearParameters ( unsigned int index) const

Returns a curvilinear representation of the parameters at 'index'.

Note
This is only available in Athena.

Definition at line 625 of file TrackParticle_v1.cxx.

625 {
626
627 static const Accessor< std::vector<float> > acc( "trackParameterCovarianceMatrices" );
628 unsigned int offset = index*15;
629 // copy the correct values into the temp matrix
630 ParametersCovMatrix_t cov;
631 auto it = acc(*this).begin()+offset;
632 Amg::expand(it,it+15,cov);
633 // retrieve the parameters to build the curvilinear frame
634 Amg::Vector3D pos(parameterX(index),parameterY(index),parameterZ(index));
636 Trk::CurvilinearParameters param(pos,mom,charge(),std::move(cov));
637
638 return param;
639 }
float parameterPX(unsigned int index) const
Returns the parameter x momentum component, for 'index'.
float parameterX(unsigned int index) const
Returns the parameter x position, for 'index'.
float parameterPY(unsigned int index) const
Returns the parameter y momentum component, for 'index'.
float parameterZ(unsigned int index) const
Returns the parameter z position, for 'index'.
float parameterY(unsigned int index) const
Returns the parameter y position, for 'index'.
float charge() const
Returns the charge.
float parameterPZ(unsigned int index) const
Returns the parameter z momentum component, for 'index'.
void expand(std::vector< float >::const_iterator it, std::vector< float >::const_iterator, AmgSymMatrix(N) &covMatrix)
Eigen::Matrix< double, 3, 1 > Vector3D
str index
Definition DeMoScan.py:362
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
static const SG::AuxElement::Accessor< ElementLink< IParticleContainer > > acc("originalObjectLink")
Object used for setting/getting the dynamic decoration in question.

◆ d0()

float xAOD::TrackParticle_v1::d0 ( ) const

Returns the \(d_0\) parameter.

◆ definingParameters()

DefiningParameters_t xAOD::TrackParticle_v1::definingParameters ( ) const

Returns a SVector of the Perigee track parameters.

i.e. a vector of \(\left(\begin{array}{c}d_0\\z_0\\\phi_0\\\theta\\q/p\end{array}\right)\)

Definition at line 175 of file TrackParticle_v1.cxx.

175 {
176 DefiningParameters_t tmp;
177 tmp << d0() , z0() , phi0() , theta() , qOverP();
178 return tmp;
179 }
float z0() const
Returns the parameter.
float theta() const
Returns the parameter, which has range 0 to .
float d0() const
Returns the parameter.
float phi0() const
Returns the parameter, which has range to .

◆ definingParametersCovMatrix()

const xAOD::ParametersCovMatrix_t xAOD::TrackParticle_v1::definingParametersCovMatrix ( ) const

Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.

Definition at line 260 of file TrackParticle_v1.cxx.

260 {
261
262 std::span<const float> covMatrixDiag;
263 if( accCovMatrixDiag.isAvailable( *this ))
264 covMatrixDiag = std::span<const float>( accCovMatrixDiag( *this ) );
265 std::span<const float> covMatrixOffDiag;
266 if( accCovMatrixOffDiag.isAvailable( *this ))
267 covMatrixOffDiag = std::span<const float>( accCovMatrixOffDiag( *this ) );
268 bool valid = true;
269 auto result = TrackingDetails::definingParametersCovMatrix( covMatrixDiag, covMatrixOffDiag, valid );
271 return result;
272 }
list valid
Definition calibdata.py:44
void covarianceUnsetHook()
Function that would be possible to use to debug what client is trying to access covariance matrix fro...
xAOD::ParametersCovMatrix_t definingParametersCovMatrix(std::span< const float > covMatrixDiag, std::span< const float > covMatrixOffDiag, bool &valid)
static const SG::AuxElement::Accessor< std::vector< float > > accCovMatrixDiag("definingParametersCovMatrixDiag")

◆ definingParametersCovMatrixDiagVec()

const std::vector< float > & xAOD::TrackParticle_v1::definingParametersCovMatrixDiagVec ( ) const

Returns the diagonal elements of the defining parameters covariance matrix.

Definition at line 326 of file TrackParticle_v1.cxx.

326 {
327
328 return accCovMatrixDiag( *this );
329 }

◆ definingParametersCovMatrixFilled()

ParametersCovMatrixFilled_t xAOD::TrackParticle_v1::definingParametersCovMatrixFilled ( ) const

Returns a 5x5 matrix describing which elements of the covariance matrix are known.

Definition at line 274 of file TrackParticle_v1.cxx.

274 {
275
276 // Create the result matrix.
278 result.setZero();
279
280 // Check if the diagonal elements are available.
281 if( accCovMatrixDiag.isAvailable( *this ) &&
282 ( static_cast< int >( accCovMatrixDiag( *this ).size() ) == result.rows() ) ) {
283
284 result.setIdentity();
285 }
286
287 bool offDiagCompr = definingParametersCovMatrixOffDiagCompr();
288
289 if(!offDiagCompr){
290
291 // Check if the off-diagonal elements are available.
292 if( accCovMatrixOffDiag.isAvailable( *this ) &&
293 ( static_cast< int >( accCovMatrixOffDiag( *this ).size() ) ==
294 ( ( result.rows() * ( result.rows() - 1 ) ) / 2 ) ) ) {
295
296 for( int i = 1; i < result.rows(); ++i ) {
297 for( int j = 0; j < i; ++j ) {
298 result.fillSymmetric( i, j, true );
299 }
300 }
301 }
302
303 }
304
305 else{
306
307 if( accCovMatrixOffDiag.isAvailable( *this ) &&
308 ( static_cast< int >( accCovMatrixOffDiag( *this ).size() ) == COVMATRIX_OFFDIAG_VEC_COMPR_SIZE ) ){
309
310 const covMatrixIndexPairVec& vecPairIndex = covMatrixComprIndexPairs();
311
312 for(const auto& pairIndex : vecPairIndex){
313 covMatrixIndex i = pairIndex.first;
314 covMatrixIndex j = pairIndex.second;
315 result.fillSymmetric( i, j, true );
316 }
317
318 }
319
320 }
321
322 // Return the object.
323 return result;
324 }
bool definingParametersCovMatrixOffDiagCompr() const
Eigen::Matrix< bool, 5, 5, 0, 5, 5 > ParametersCovMatrixFilled_t

◆ definingParametersCovMatrixOffDiagCompr()

bool xAOD::TrackParticle_v1::definingParametersCovMatrixOffDiagCompr ( ) const

Definition at line 380 of file TrackParticle_v1.cxx.

380 {
381
382 std::span<const float> covMatrixOffDiag;
383 if( accCovMatrixOffDiag.isAvailable( *this ))
384 covMatrixOffDiag = std::span<const float>( accCovMatrixOffDiag( *this ) );
386 }
bool definingParametersCovMatrixOffDiagCompr(std::span< const float > covMatrixOffDiag)

◆ definingParametersCovMatrixOffDiagVec()

const std::vector< float > & xAOD::TrackParticle_v1::definingParametersCovMatrixOffDiagVec ( ) const

Returns the correlation coefficient associated with the off-diagonal elements of the covariance matrix = cov(X,Y)/sqrt(cov(X,X)*cov(Y,Y)).

Definition at line 331 of file TrackParticle_v1.cxx.

331 {
332
333 return accCovMatrixOffDiag( *this );
334 }

◆ definingParametersCovMatrixVec()

std::vector< float > xAOD::TrackParticle_v1::definingParametersCovMatrixVec ( ) const

Returns the length 6 vector containing the elements of defining parameters covariance matrix.

Definition at line 336 of file TrackParticle_v1.cxx.

336 {
337
338 std::vector< float > vec;
340 Amg::compress(cov,vec);
341 return vec;
342
343 }
std::vector< size_t > vec
#define AmgSymMatrix(dim)

◆ e()

double xAOD::TrackParticle_v1::e ( ) const
finaloverridevirtual

The total energy of the particle.

Implements xAOD::IParticle.

Definition at line 111 of file TrackParticle_v1.cxx.

111 {
112 return genvecP4().E();
113 }
GenVecFourMom_t genvecP4() const
The full 4-momentum of the particle : GenVector form.

◆ eta()

double xAOD::TrackParticle_v1::eta ( ) const
finaloverridevirtual

The pseudorapidity ( \(\eta\)) of the particle.

Implements xAOD::IParticle.

Definition at line 79 of file TrackParticle_v1.cxx.

79 {
80 return genvecP4().Eta();
81 }

◆ genvecP4()

TrackParticle_v1::GenVecFourMom_t xAOD::TrackParticle_v1::genvecP4 ( ) const

The full 4-momentum of the particle : GenVector form.

Definition at line 118 of file TrackParticle_v1.cxx.

118 {
119 return TrackingDetails::genvecP4(qOverP(), theta(), phi(), m());
120 }
virtual double m() const override final
The invariant mass of the particle..
GenVecFourMom_t genvecP4(float qOverP, float thetaT, float phiT, double m)

◆ hasValidTime()

uint8_t xAOD::TrackParticle_v1::hasValidTime ( ) const

Returns whether or not the track has a valid time.

◆ hitPattern()

uint32_t xAOD::TrackParticle_v1::hitPattern ( ) const

◆ identifierOfFirstHit()

uint64_t xAOD::TrackParticle_v1::identifierOfFirstHit ( ) const

Returns the offline identifier of the first hit.

◆ indexOfParameterAtPosition()

bool xAOD::TrackParticle_v1::indexOfParameterAtPosition ( unsigned int & index,
ParameterPosition position ) const

Function to determine if this TrackParticle contains track parameters at a certain position, and if so, what the 'index' is.

Parameters
[in]indexFilled with the index of the track parameters at 'position' - untouched otherwise.
[out]positionThe location in the detector of the required track parameters.
Returns
Returns 'true' if the TrackParticle parameters at 'position', returns False otherwise.

Definition at line 605 of file TrackParticle_v1.cxx.

606 {
607 size_t maxParameters = numberOfParameters();
608 bool foundParameters=false;
609 for (size_t i=0; i<maxParameters; ++i){
610 if (parameterPosition(i)==position){
611 foundParameters=true;
612 index=i;
613 break;
614 }
615 }
616 return foundParameters;
617 }
xAOD::ParameterPosition parameterPosition(unsigned int index) const
Return the ParameterPosition of the parameters at 'index'.
size_t numberOfParameters() const
Returns the number of additional parameters stored in the TrackParticle.

◆ isAvailable()

template<class T>
XAOD_AUXDATA_DEPRECATED bool xAOD::IParticle::isAvailable ( const std::string & name,
const std::string & clsname = "" ) const
inlineinherited

Check if a user property is available for reading or not.

This function should be used to check if a user property which may or may not exist, is set on the object.

Parameters
nameName of the auxiliary variable
clsnameThe name of the associated class. May be blank
Returns
Whether the decoration exists or not

Definition at line 135 of file Event/xAOD/xAODBase/xAODBase/IParticle.h.

136 {
137
138 return SG::ConstAccessor< T >(name, clsname).isAvailable(*this);
139 }
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.

◆ isAvailableWritable()

template<class T>
XAOD_AUXDATA_DEPRECATED bool xAOD::IParticle::isAvailableWritable ( const std::string & name,
const std::string & clsname = "" ) const
inlineinherited

Check if a user property is available for writing or not.

This function can be used to check whether it will be possible to set a user property on the object.

Parameters
nameName of the auxiliary variable
clsnameThe name of the associated class. May be blank
Returns
Whether the decoration is possible to set

Definition at line 152 of file Event/xAOD/xAODBase/xAODBase/IParticle.h.

153 {
154
155 return SG::Accessor< T >(name, clsname).isAvailableWritable(*this);
156 }
bool isAvailableWritable(ELT &e) const
Test to see if this variable exists in the store and is writable.

◆ m()

double xAOD::TrackParticle_v1::m ( ) const
finaloverridevirtual

The invariant mass of the particle..

Implements xAOD::IParticle.

Definition at line 85 of file TrackParticle_v1.cxx.

85 {
86 // Codes using a fitter set a hypothesis, and the
87 // particular fitter that was employed..
88 // A mass is never set/stored.
89 //
90 // In the past we were returning the mass of a charged pion always
91 //
92 // This created a confusion on why TrackParticles created by
93 // specific lepton fitter have a pion mass (the leptons per se have the
94 // correct mass). Lets try to remedy this.
96 if (hypo == xAOD::electron) {
97 // Since GX2 also set sometimes the hypo to electron
98 // lets also check for GSF.
100 if (fitter == xAOD::GaussianSumFilter) {
102 }
103 }
104 if (hypo == xAOD::muon) {
106 }
107 // default charged pion
109 }
TrackFitter trackFitter() const
Returns the fitter.
ParticleHypothesis particleHypothesis() const
Returns the particle hypothesis used for Track fitting.
const ShapeFitter * fitter
constexpr double muonMassInMeV
the mass of the muon (in MeV)
constexpr double chargedPionMassInMeV
the mass of the charged pion (in MeV)
constexpr double electronMassInMeV
the mass of the electron (in MeV)
@ GaussianSumFilter
Tracks from Gaussian Sum Filter.

◆ numberDoF()

float xAOD::TrackParticle_v1::numberDoF ( ) const

Returns the number of degrees of freedom of the overall track or vertex fit as float.

◆ numberOfIBLOverflowsdEdx()

uint8_t xAOD::TrackParticle_v1::numberOfIBLOverflowsdEdx ( ) const

◆ numberOfParameters()

size_t xAOD::TrackParticle_v1::numberOfParameters ( ) const

Returns the number of additional parameters stored in the TrackParticle.

number of parameters should be the size of positions we need for them

Definition at line 505 of file TrackParticle_v1.cxx.

505 {
507 static const Accessor< std::vector<uint8_t> > acc( "parameterPosition" );
508 if(! acc.isAvailable( *this )) return 0;
509 return acc(*this).size();
510 }

◆ numberOfUsedHitsdEdx()

uint8_t xAOD::TrackParticle_v1::numberOfUsedHitsdEdx ( ) const

◆ operator=()

TrackParticle_v1 & xAOD::TrackParticle_v1::operator= ( const TrackParticle_v1 & tp)

Assignment operator. This can involve creating and copying an Auxilary store, and so should be used sparingly.

Definition at line 59 of file TrackParticle_v1.cxx.

59 {
60 if(this == &tp) return *this;
61
62 if( ( ! hasStore() ) && ( ! container() ) ) {
64 }
65 this->IParticle::operator=( tp );
66#ifndef XAOD_ANALYSIS
67 // assume that this copy will create new cache as needed
68 m_perigeeParameters.reset();
69#endif // not XAOD_ANALYSIS
70 return *this;
71 }
IParticle & operator=(const IParticle &)=default
CxxUtils::CachedValue< Trk::Perigee > m_perigeeParameters
Cached MeasuredPerigee, built from this object.
const SG::AuxVectorData * container() const
Return the container holding this element.
bool hasStore() const
Return true if this object has an associated store.

◆ p4()

TrackParticle_v1::FourMom_t xAOD::TrackParticle_v1::p4 ( ) const
finaloverridevirtual

The full 4-momentum of the particle.

Implements xAOD::IParticle.

Definition at line 122 of file TrackParticle_v1.cxx.

122 {
124 using namespace std;
125 float p = 10.e6; // 10 TeV (default value for very high pt muons, with qOverP==0)
126 if (fabs(qOverP())>0.) p = 1/fabs(qOverP());
127 float thetaT = theta();
128 float phiT = phi();
129 float sinTheta= sin(thetaT);
130 float px = p*sinTheta*cos(phiT);
131 float py = p*sinTheta*sin(phiT);
132 float pz = p*cos(thetaT);
133 float e = pow (m(),2) +
134 pow( px,2) + pow( py,2) + pow( pz,2);
135 p4.SetPxPyPzE( px, py, pz, sqrt(e) );
136 return p4;
137 }
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
IParticle::FourMom_t FourMom_t
Definition of the 4-momentum type.
virtual double e() const override final
The total energy of the particle.
constexpr int pow(int x)
Definition conifer.h:27

◆ parameterPosition()

xAOD::ParameterPosition xAOD::TrackParticle_v1::parameterPosition ( unsigned int index) const

Return the ParameterPosition of the parameters at 'index'.

Definition at line 599 of file TrackParticle_v1.cxx.

600 {
601 static const Accessor< std::vector<uint8_t> > acc( "parameterPosition" );
602 return static_cast<xAOD::ParameterPosition>(acc(*this).at(index));
603 }
ParameterPosition
Enum allowing us to know where in ATLAS the parameters are defined.

◆ parameterPX()

float xAOD::TrackParticle_v1::parameterPX ( unsigned int index) const

Returns the parameter x momentum component, for 'index'.

Definition at line 564 of file TrackParticle_v1.cxx.

564 {
565 static const Accessor< std::vector<float> > acc( "parameterPX" );
566 return acc(*this).at(index);
567 }

◆ parameterPY()

float xAOD::TrackParticle_v1::parameterPY ( unsigned int index) const

Returns the parameter y momentum component, for 'index'.

Definition at line 569 of file TrackParticle_v1.cxx.

569 {
570 static const Accessor< std::vector<float> > acc( "parameterPY" );
571 return acc(*this).at(index);
572 }

◆ parameterPZ()

float xAOD::TrackParticle_v1::parameterPZ ( unsigned int index) const

Returns the parameter z momentum component, for 'index'.

Definition at line 574 of file TrackParticle_v1.cxx.

574 {
575 static const Accessor< std::vector<float> > acc( "parameterPZ" );
576 return acc(*this).at(index);
577 }

◆ parameterX()

float xAOD::TrackParticle_v1::parameterX ( unsigned int index) const

Returns the parameter x position, for 'index'.

Definition at line 549 of file TrackParticle_v1.cxx.

549 {
550 static const Accessor< std::vector<float> > acc( "parameterX" );
551 return acc(*this).at(index);
552 }

◆ parameterY()

float xAOD::TrackParticle_v1::parameterY ( unsigned int index) const

Returns the parameter y position, for 'index'.

Definition at line 554 of file TrackParticle_v1.cxx.

554 {
555 static const Accessor< std::vector<float> > acc( "parameterY" );
556 return acc(*this).at(index);
557 }

◆ parameterZ()

float xAOD::TrackParticle_v1::parameterZ ( unsigned int index) const

Returns the parameter z position, for 'index'.

Definition at line 559 of file TrackParticle_v1.cxx.

559 {
560 static const Accessor< std::vector<float> > acc( "parameterZ" );
561 return acc(*this).at(index);
562 }

◆ particleHypothesis()

xAOD::ParticleHypothesis xAOD::TrackParticle_v1::particleHypothesis ( ) const

Returns the particle hypothesis used for Track fitting.

Definition at line 680 of file TrackParticle_v1.cxx.

680 {
681 static const Accessor<uint8_t> acc("particleHypothesis");
682 if (!acc.isAvailable(*this)) {
683 return xAOD::pion;
684 }
685 return static_cast<xAOD::ParticleHypothesis>(acc(*this));
686 }

◆ patternRecoInfo()

std::bitset< xAOD::NumberOfTrackRecoInfo > xAOD::TrackParticle_v1::patternRecoInfo ( ) const

Access method for pattern recognition algorithm.

Definition at line 658 of file TrackParticle_v1.cxx.

659 {
660 static const Accessor< uint64_t > acc( "patternRecoInfo" );
661 std::bitset<xAOD::NumberOfTrackRecoInfo> tmp(acc(*this));
662 return tmp;
663 }

◆ perigeeParameters()

const Trk::Perigee & xAOD::TrackParticle_v1::perigeeParameters ( ) const

Returns the Trk::MeasuredPerigee track parameters.

These are defined as: \(\left(\begin{array}{c}d_0\\z_0\\\phi_0\\\theta\\q/p\\\end{array}\right)\)

Note
This is only available in Athena.

Definition at line 437 of file TrackParticle_v1.cxx.

437 {
438
439 // Require the cache to be valid and check if the cached pointer has been set
440 if(m_perigeeParameters.isValid()){
441 return *(m_perigeeParameters.ptr());
442 }
443 static const Accessor< float > acc1( "d0" );
444 static const Accessor< float > acc2( "z0" );
445 static const Accessor< float > acc3( "phi" );
446 static const Accessor< float > acc4( "theta" );
447 static const Accessor< float > acc5( "qOverP" );
448 static const Accessor< std::vector<float> > acc6( "definingParametersCovMatrix" );
449 ParametersCovMatrix_t cov = ParametersCovMatrix_t(definingParametersCovMatrix());
450 static const Accessor< float > acc7( "beamlineTiltX" );
451 static const Accessor< float > acc8( "beamlineTiltY" );
452
453 if(!acc7.isAvailable( *this ) || !acc8.isAvailable( *this )){
454 Trk::Perigee tmpPerigeeParameters(
455 acc1(*this),
456 acc2(*this),
457 acc3(*this),
458 acc4(*this),
459 acc5(*this),
460 Trk::PerigeeSurface(Amg::Vector3D(vx(), vy(), vz())),
461 std::move(cov));
462 m_perigeeParameters.set(tmpPerigeeParameters);
463 return *(m_perigeeParameters.ptr());
464 }
465
466 Amg::Translation3D amgtranslation(vx(),vy(),vz());
467 Amg::Transform3D pAmgTransf = amgtranslation * Amg::RotationMatrix3D::Identity();
468 pAmgTransf *= Amg::AngleAxis3D(acc8(*this), Amg::Vector3D(0.,1.,0.));
469 pAmgTransf *= Amg::AngleAxis3D(acc7(*this), Amg::Vector3D(1.,0.,0.));
470 Trk::Perigee tmpPerigeeParameters(acc1(*this),
471 acc2(*this),
472 acc3(*this),
473 acc4(*this),
474 acc5(*this),
475 pAmgTransf,
476 std::move(cov));
477
478 m_perigeeParameters.set(tmpPerigeeParameters);
479 return *(m_perigeeParameters.ptr());
480 }
float vx() const
The x origin for the parameters.
float vz() const
The z origin for the parameters.
float vy() const
The y origin for the parameters.
Eigen::AngleAxisd AngleAxis3D
Eigen::Affine3d Transform3D
Eigen::Translation< double, 3 > Translation3D
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee

◆ phi()

virtual double xAOD::TrackParticle_v1::phi ( ) const
finaloverridevirtual

The azimuthal angle ( \(\phi\)) of the particle (has range \(-\pi\) to \(+\pi\).).

Implements xAOD::IParticle.

◆ phi0()

float xAOD::TrackParticle_v1::phi0 ( ) const

Returns the \(\phi\) parameter, which has range \(-\pi\) to \(+\pi\).

Definition at line 150 of file TrackParticle_v1.cxx.

150 {
151
152 static const Accessor< float > acc( "phi" );
153 return acc( *this );
154 }

◆ pt()

double xAOD::TrackParticle_v1::pt ( ) const
finaloverridevirtual

The transverse momentum ( \(p_T\)) of the particle.

Implements xAOD::IParticle.

Definition at line 75 of file TrackParticle_v1.cxx.

75 {
76 return genvecP4().Pt();
77 }

◆ qOverP()

float xAOD::TrackParticle_v1::qOverP ( ) const

Returns the \(q/p\) parameter.

◆ radiusOfFirstHit()

float xAOD::TrackParticle_v1::radiusOfFirstHit ( ) const

Returns the radius of the first hit.

◆ rapidity()

double xAOD::TrackParticle_v1::rapidity ( ) const
finaloverridevirtual

The true rapidity (y) of the particle.

Implements xAOD::IParticle.

Definition at line 114 of file TrackParticle_v1.cxx.

114 {
115 return genvecP4().Rapidity();
116 }

◆ resetCache()

void xAOD::TrackParticle_v1::resetCache ( )

Reset the internal cache of the object.

Definition at line 772 of file TrackParticle_v1.cxx.

772 {
773#ifndef XAOD_ANALYSIS
774 m_perigeeParameters.reset();
775#endif // not XAOD_ANALYSIS
776 }

◆ setBeamlineTiltX()

void xAOD::TrackParticle_v1::setBeamlineTiltX ( float tiltX)

◆ setBeamlineTiltY()

void xAOD::TrackParticle_v1::setBeamlineTiltY ( float tiltY)

◆ setDefiningParameters() [1/2]

void xAOD::TrackParticle_v1::setDefiningParameters ( float d0,
float z0,
float phi0,
float theta,
float qOverP )

Set the defining parameters.

Definition at line 181 of file TrackParticle_v1.cxx.

181 {
182#ifndef XAOD_ANALYSIS
183 // reset perigee cache if existing
184 if(m_perigeeParameters.isValid()) {
185 m_perigeeParameters.reset();
186 }
187#endif // not XAOD_ANALYSIS
188 static const Accessor< float > acc1( "d0" );
189 acc1( *this ) = d0;
190
191 static const Accessor< float > acc2( "z0" );
192 acc2( *this ) = z0;
193
194 static const Accessor< float > acc3( "phi" );
195 acc3( *this ) = phi0;
196
197 static const Accessor< float > acc4( "theta" );
198 acc4( *this ) = theta;
199
200 static const Accessor< float > acc5( "qOverP" );
201 acc5( *this ) = qOverP;
202
203 return;
204 }

◆ setDefiningParameters() [2/2]

void xAOD::TrackParticle_v1::setDefiningParameters ( float d0,
float z0,
float phi0,
float theta,
float qOverP,
float time )

Definition at line 206 of file TrackParticle_v1.cxx.

206 {
208 setTime(time);
209 return;
210 }
float time() const
Returns the time.
void setDefiningParameters(float d0, float z0, float phi0, float theta, float qOverP)
Set the defining parameters.

◆ setDefiningParametersCovMatrix()

void xAOD::TrackParticle_v1::setDefiningParametersCovMatrix ( const ParametersCovMatrix_t & cov)

Set the defining parameters covariance matrix.

Definition at line 227 of file TrackParticle_v1.cxx.

227 {
228
229#ifndef XAOD_ANALYSIS
230 // reset perigee cache if existing
231 if(m_perigeeParameters.isValid()) {
232 m_perigeeParameters.reset();
233 }
234#endif // not XAOD_ANALYSIS
235
236 // Extract the diagonal elements from the matrix.
237 std::vector< float > diagVec;
238 diagVec.reserve( cov.rows() );
239 for( int i = 0; i < cov.rows(); ++i ) {
240 diagVec.push_back( cov( i, i ) );
241 }
242 // Set the variable.
244
245 // Extract the off-diagonal elements from the matrix.
246 std::vector< float > offDiagVec;
247 offDiagVec.reserve( ( ( cov.rows() - 1 ) * cov.rows() ) / 2 );
248 for( int i = 1; i < cov.rows(); ++i ) {
249 for( int j = 0; j < i; ++j ) {
250 float offDiagCoeff = (cov( i, i )>0 && cov( j, j )>0) ? cov( i, j )/sqrt(cov( i, i )*cov( j, j )) : 0;
251 offDiagVec.push_back( offDiagCoeff );
252 }
253 }
254 // Set the variable.
256
257 return;
258 }
void setDefiningParametersCovMatrixOffDiagVec(const std::vector< float > &vec)
Set the off-diagonal elements of the defining parameters covariance matrix.
void setDefiningParametersCovMatrixDiagVec(const std::vector< float > &vec)
Set the defining parameters covariance matrix using a length 15 vector.

◆ setDefiningParametersCovMatrixDiagVec()

void xAOD::TrackParticle_v1::setDefiningParametersCovMatrixDiagVec ( const std::vector< float > & vec)

Set the defining parameters covariance matrix using a length 15 vector.

Set the diagonal elements of the defining parameters covariance matrix

Definition at line 345 of file TrackParticle_v1.cxx.

345 {
346
347 if (vec.size() != ParametersCovMatrix_t::RowsAtCompileTime) {
348 throw std::runtime_error(
349 "Setting track definingParametersCovMatrixDiag with vector of size " +
350 std::to_string(vec.size()) + " instead of expected " +
351 std::to_string(ParametersCovMatrix_t::RowsAtCompileTime) +
352 " is not supported");
353 }
354
355 accCovMatrixDiag( *this ) = vec;
356 return;
357 }

◆ setDefiningParametersCovMatrixOffDiagVec()

void xAOD::TrackParticle_v1::setDefiningParametersCovMatrixOffDiagVec ( const std::vector< float > & vec)

Set the off-diagonal elements of the defining parameters covariance matrix.

Definition at line 359 of file TrackParticle_v1.cxx.

359 {
360
361 bool offDiagCompr = definingParametersCovMatrixOffDiagCompr();
362
363 unsigned int uncompr_size = ( ( ( ParametersCovMatrix_t::RowsAtCompileTime - 1 ) *
364 ParametersCovMatrix_t::RowsAtCompileTime ) / 2 );
365 unsigned int size = offDiagCompr ? COVMATRIX_OFFDIAG_VEC_COMPR_SIZE : uncompr_size;
366
367 if( !(vec.size() == size || vec.size() == uncompr_size) ){ //If off-diagonal elements are already compressed, can either set with uncompressed or compressed vector
368 throw std::runtime_error(
369 "Setting track definingParametersCovMatrixOffDiag with vector of "
370 "size " +
371 std::to_string(vec.size()) + " instead of expected " +
372 std::to_string(size) + " or " + std::to_string(uncompr_size) +
373 " is not supported");
374 }
375
376 accCovMatrixOffDiag( *this ) = vec;
377 return;
378 }
size_t size() const
Number of registered mappings.

◆ setDefiningParametersCovMatrixVec()

void xAOD::TrackParticle_v1::setDefiningParametersCovMatrixVec ( const std::vector< float > & cov)

Definition at line 412 of file TrackParticle_v1.cxx.

412 {
413
414 xAOD::ParametersCovMatrix_t covMatrix;
415 if( !cov.empty() ) Amg::expand( cov.begin(), cov.end(),covMatrix );
416 else covMatrix.setIdentity();
418
419 }
void setDefiningParametersCovMatrix(const ParametersCovMatrix_t &cov)
Set the defining parameters covariance matrix.

◆ setFitQuality()

void xAOD::TrackParticle_v1::setFitQuality ( float chiSquared,
float numberDoF )

Set the 'Fit Quality' information.

Definition at line 486 of file TrackParticle_v1.cxx.

486 {
487 static const Accessor< float > acc1( "chiSquared" );
488 acc1( *this ) = chiSquared;
489 static const Accessor< float > acc2( "numberDoF" );
490 acc2( *this ) = numberDoF;
491 }
float numberDoF() const
Returns the number of degrees of freedom of the overall track or vertex fit as float.
float chiSquared() const
Returns the of the overall track fit.

◆ setHasValidTime()

void xAOD::TrackParticle_v1::setHasValidTime ( uint8_t hasValidTime)

◆ setHitPattern()

void xAOD::TrackParticle_v1::setHitPattern ( uint32_t hitpattern)

◆ setIdentifierOfFirstHit()

void xAOD::TrackParticle_v1::setIdentifierOfFirstHit ( uint64_t id)

Set the offline identifier of the first hit.

◆ setNumberOfIBLOverflowsdEdx()

void xAOD::TrackParticle_v1::setNumberOfIBLOverflowsdEdx ( uint8_t numoverflows)

◆ setNumberOfUsedHitsdEdx()

void xAOD::TrackParticle_v1::setNumberOfUsedHitsdEdx ( uint8_t numhits)

◆ setParameterPosition()

void xAOD::TrackParticle_v1::setParameterPosition ( unsigned int index,
xAOD::ParameterPosition pos )

Set the 'position' (i.e. where it is in ATLAS) of the parameter at 'index', using the ParameterPosition enum.

Definition at line 619 of file TrackParticle_v1.cxx.

619 {
620 static const Accessor< std::vector<uint8_t> > acc( "parameterPosition" );
621 acc( *this ).at(index) = static_cast<uint8_t>(pos);
622 }

◆ setParametersOrigin()

void xAOD::TrackParticle_v1::setParametersOrigin ( float x,
float y,
float z )

Set the origin for the parameters.

Definition at line 425 of file TrackParticle_v1.cxx.

425 {
426 static const Accessor< float > acc1( "vx" );
427 acc1( *this ) = x;
428
429 static const Accessor< float > acc2( "vy" );
430 acc2( *this ) = y;
431
432 static const Accessor< float > acc3( "vz" );
433 acc3( *this ) = z;
434 }
#define y
#define x
#define z

◆ setParticleHypothesis()

void xAOD::TrackParticle_v1::setParticleHypothesis ( const ParticleHypothesis hypo)

Method for setting the particle type, using the ParticleHypothesis enum.

Definition at line 675 of file TrackParticle_v1.cxx.

675 {
676 static const Accessor<uint8_t> acc("particleHypothesis");
677 acc(*this) = static_cast<uint8_t>(value);
678 }

◆ setPatternRecognitionInfo() [1/2]

void xAOD::TrackParticle_v1::setPatternRecognitionInfo ( const std::bitset< xAOD::NumberOfTrackRecoInfo > & patternReco)

Method setting the pattern recognition algorithm, using a bitset.

The bitset should be created using the TrackPatternRecoInfo enum as follows:

const std::bitset<xAOD::NumberOfTrackRecoInfo> patternReco;
patternReco.set(xAOD::Fatras);
@ Fatras
Track from FATRAS.

Definition at line 670 of file TrackParticle_v1.cxx.

670 {
671 static const Accessor< uint64_t > acc( "patternRecoInfo" );
672 acc( *this ) = patternReco.to_ullong();
673 }

◆ setPatternRecognitionInfo() [2/2]

void xAOD::TrackParticle_v1::setPatternRecognitionInfo ( uint64_t patternReco)

Method setting the pattern recognition algorithm, using a 64-bit int (which is faster than using a bitset).

The bit set should be created using the TrackPatternRecoInfo enum as follows:

uint64_t patternReco;
patternReco |= (1<<static_cast<uint64_t>(xAOD::Fatras))

Definition at line 665 of file TrackParticle_v1.cxx.

665 {
666 static const Accessor< uint64_t > acc( "patternRecoInfo" );
667 acc( *this ) = patternReco;
668 }

◆ setRadiusOfFirstHit()

void xAOD::TrackParticle_v1::setRadiusOfFirstHit ( float radius)

Set the radius of the first hit.

◆ setSummaryValue() [1/2]

void xAOD::TrackParticle_v1::setSummaryValue ( float & value,
const SummaryType & information )

Set method for TrackSummary values.

Definition at line 710 of file TrackParticle_v1.cxx.

710 {
711 const xAOD::TrackParticle_v1::Accessor< float >* acc = trackSummaryAccessorV1<float>( information );
712 // Set the value:
713 ( *acc )( *this ) = value;
714 }
const SG::AuxElement::Accessor< float > * trackSummaryAccessorV1< float >(xAOD::SummaryType type)

◆ setSummaryValue() [2/2]

void xAOD::TrackParticle_v1::setSummaryValue ( uint8_t & value,
const SummaryType & information )

Set method for TrackSummary values.

Definition at line 704 of file TrackParticle_v1.cxx.

704 {
705 const xAOD::TrackParticle_v1::Accessor< uint8_t >* acc = trackSummaryAccessorV1<uint8_t>( information );
706 // Set the value:
707 ( *acc )( *this ) = value;
708 }
const SG::AuxElement::Accessor< uint8_t > * trackSummaryAccessorV1< uint8_t >(xAOD::SummaryType type)

◆ setTime()

void xAOD::TrackParticle_v1::setTime ( float time)

Definition at line 212 of file TrackParticle_v1.cxx.

212 {
213 static const SG::AuxElement::Accessor< float > acc("time");
214 acc( *this ) = time;
215 }

◆ setTimeResolution()

void xAOD::TrackParticle_v1::setTimeResolution ( float timeResolution)

Definition at line 217 of file TrackParticle_v1.cxx.

217 {
218 static const SG::AuxElement::Accessor< float > acc("timeResolution");
219 acc( *this ) = timeRes;
220 }

◆ setTrackFitter()

void xAOD::TrackParticle_v1::setTrackFitter ( const TrackFitter fitter)

Method for setting the fitter, using the TrackFitter enum.

Definition at line 645 of file TrackParticle_v1.cxx.

645 {
646 static const Accessor<uint8_t> acc("trackFitter");
647 acc(*this) = static_cast<uint8_t>(value);
648 }

◆ setTrackLink()

void xAOD::TrackParticle_v1::setTrackLink ( const ElementLink< TrackCollection > & track)

Set the link to the original track.

Note
This is only available in Athena.

Definition at line 745 of file TrackParticle_v1.cxx.

746 {
747
748 // The accessor:
749 static const Accessor< ElementLink< TrackCollection > > acc( "trackLink" );
750
751 // Do the deed:
752 acc( *this ) = el;
753 return;
754 }

◆ setTrackParameterCovarianceMatrix()

void xAOD::TrackParticle_v1::setTrackParameterCovarianceMatrix ( unsigned int index,
std::vector< float > & cov )

Set the cov matrix of the parameter at 'index', using a vector of floats.

The vector \(\mathrm{v}(a1,a2,a3 ... a_{15})\) represents the lower diagonal, i.e. it gives a matrix of \(\left(\begin{array}{ccccc} a_1 & a_2 & a_4 & a_7 & a_{11} \\ a_2 & a_3 & a_5 & a_8 & a_{12} \\ a_4 & a_5 & a_6 & a_9 & a_{13} \\ a_7 & a_8 & a_9 & a_{10} & a_{14} \\ a_{11} & a_{12} & a_{13} & a_{14} & a_{15} \end{array}\right)\)

Definition at line 590 of file TrackParticle_v1.cxx.

590 {
591 assert(cov.size()==15);
592 unsigned int offset = index*15;
593 static const Accessor< std::vector < float > > acc( "trackParameterCovarianceMatrices" );
594 std::vector<float>& v = acc(*this);
595 v.resize(offset+15);
596 std::copy(cov.begin(),cov.end(),v.begin()+offset );
597 }

◆ setTrackParameters()

void xAOD::TrackParticle_v1::setTrackParameters ( std::vector< std::vector< float > > & parameters)

Set the parameters via the passed vector of vectors.

The vector<float> should be of size 6: x,y,z,px,py,pz (charge is stored elsewhere)

Definition at line 519 of file TrackParticle_v1.cxx.

519 {
520 static const Accessor< std::vector < float > > acc1( "parameterX" );
521 static const Accessor< std::vector < float > > acc2( "parameterY" );
522 static const Accessor< std::vector < float > > acc3( "parameterZ" );
523 static const Accessor< std::vector < float > > acc4( "parameterPX" );
524 static const Accessor< std::vector < float > > acc5( "parameterPY" );
525 static const Accessor< std::vector < float > > acc6( "parameterPZ" );
526 static const Accessor< std::vector<uint8_t> > acc7( "parameterPosition" );
527
528 acc1(*this).resize(parameters.size());
529 acc2(*this).resize(parameters.size());
530 acc3(*this).resize(parameters.size());
531 acc4(*this).resize(parameters.size());
532 acc5(*this).resize(parameters.size());
533 acc6(*this).resize(parameters.size());
534 acc7(*this).resize(parameters.size());
535
536 unsigned int index=0;
537 std::vector<std::vector<float> >::const_iterator it=parameters.begin(), itEnd=parameters.end();
538 for (;it!=itEnd;++it,++index){
539 assert((*it).size()==6);
540 acc1(*this).at(index)=(*it).at(0);
541 acc2(*this).at(index)=(*it).at(1);
542 acc3(*this).at(index)=(*it).at(2);
543 acc4(*this).at(index)=(*it).at(3);
544 acc5(*this).at(index)=(*it).at(4);
545 acc6(*this).at(index)=(*it).at(5);
546 }
547 }

◆ setTrackProperties()

void xAOD::TrackParticle_v1::setTrackProperties ( const TrackProperties properties)

Methods setting the TrackProperties.

◆ summaryValue() [1/2]

bool xAOD::TrackParticle_v1::summaryValue ( float & value,
const SummaryType & information ) const

Accessor for TrackSummary values.

If 'information' is stored in this TrackParticle and is of the correct templated type T, then the function fills 'value' and returns 'true', otherwise it returns 'false', and does not touch 'value'. See below for an example of how this is intended to be used.

if( myParticle.summaryValue(numberOfBLayerHits,xAOD::numberOfBLayerHits) ){
ATH_MSG_INFO("Successfully retrieved the integer value, numberOfBLayerHits");
}
float numberOfCscPhiHits=0.0; //Wrong! This is actually an int too.
if( !myParticle.summaryValue(numberOfCscPhiHits,xAOD::numberOfCscPhiHits) ){
ATH_MSG_INFO("Types must match!");
}
#define ATH_MSG_INFO(x)
@ numberOfBLayerHits
these are the hits in the first pixel layer, i.e.
Parameters
[in]informationThe information being requested. This is not guaranteed to be stored in all TrackParticles.
[out]valueOnly filled if this TrackParticle contains 'information', and the types match.
Returns
Returns 'true' if the TrackParticle contains 'information', and its concrete type matches 'value' (templated type T).

Definition at line 696 of file TrackParticle_v1.cxx.

696 {
697 const xAOD::TrackParticle_v1::Accessor< float >* acc = trackSummaryAccessorV1<float>( information );
698 if( ( ! acc ) || ( ! acc->isAvailable( *this ) ) ) return false;
699 // Retrieve the value:
700 value = ( *acc )( *this );
701 return true;
702 }

◆ summaryValue() [2/2]

bool xAOD::TrackParticle_v1::summaryValue ( uint8_t & value,
const SummaryType & information ) const

Accessor for TrackSummary values.

If 'information' is stored in this TrackParticle and is of the correct templated type T, then the function fills 'value' and returns 'true', otherwise it returns 'false', and does not touch 'value'. See below for an example of how this is intended to be used.

if( myParticle.summaryValue(numberOfBLayerHits,xAOD::numberOfBLayerHits) ){
ATH_MSG_INFO("Successfully retrieved the integer value, numberOfBLayerHits");
}
float numberOfCscPhiHits=0.0; //Wrong! This is actually an int too.
if( !myParticle.summaryValue(numberOfCscPhiHits,xAOD::numberOfCscPhiHits) ){
ATH_MSG_INFO("Types must match!");
}
Parameters
[in]informationThe information being requested. This is not guaranteed to be stored in all TrackParticles.
[out]valueOnly filled if this TrackParticle contains 'information', and the types match.
Returns
Returns 'true' if the TrackParticle contains 'information', and its concrete type matches 'value' (templated type T).

Definition at line 688 of file TrackParticle_v1.cxx.

688 {
689 const xAOD::TrackParticle_v1::Accessor< uint8_t >* acc = trackSummaryAccessorV1<uint8_t>( information );
690 if( ( ! acc ) || ( ! acc->isAvailable( *this ) ) ) return false;
691 // Retrieve the value:
692 value = ( *acc )( *this );
693 return true;
694 }

◆ theta()

float xAOD::TrackParticle_v1::theta ( ) const

Returns the \(\theta\) parameter, which has range 0 to \(\pi\).

◆ time()

float xAOD::TrackParticle_v1::time ( ) const

Returns the time.

Definition at line 161 of file TrackParticle_v1.cxx.

161 {
162 static const SG::AuxElement::Accessor< uint8_t > acc("hasValidTime");
163 if( !acc.isAvailable( *this) || !static_cast<bool>(hasValidTime()) ) throw std::runtime_error( "Unavailable TrackParticle time requested" );
164 static const SG::AuxElement::Accessor< float > accTime("time");
165 return accTime( *this );
166 }
uint8_t hasValidTime() const
Returns whether or not the track has a valid time.

◆ timeResolution()

float xAOD::TrackParticle_v1::timeResolution ( ) const

Returns the time resolution.

Definition at line 168 of file TrackParticle_v1.cxx.

168 {
169 static const SG::AuxElement::Accessor< uint8_t > acc("hasValidTime");
170 if( !acc.isAvailable( *this) || !static_cast<bool>(hasValidTime()) ) throw std::runtime_error( "Unavailable TrackParticle timeResolution requested" );
171 static const SG::AuxElement::Accessor< float > accTimeRes("timeResolution");
172 return accTimeRes( *this );
173 }

◆ track()

const Trk::Track * xAOD::TrackParticle_v1::track ( ) const

Returns a pointer (which can be NULL) to the Trk::Track which was used to make this TrackParticle.

Note
This is only available in Athena.

Definition at line 756 of file TrackParticle_v1.cxx.

756 {
757
758 // The accessor:
759 static const ConstAccessor< ElementLink< TrackCollection > > acc( "trackLink" );
760
761 if( ! acc.isAvailable( *this ) ) {
762 return nullptr;
763 }
764 if( ! acc( *this ).isValid() ) {
765 return nullptr;
766 }
767
768 return *( acc( *this ) );
769 }

◆ trackFitter()

xAOD::TrackFitter xAOD::TrackParticle_v1::trackFitter ( ) const

Returns the fitter.

Definition at line 650 of file TrackParticle_v1.cxx.

650 {
651 static const Accessor<uint8_t> acc("trackFitter");
652 if (!acc.isAvailable(*this)) {
654 }
655 return static_cast<xAOD::TrackFitter>(acc(*this));
656 }
TrackFitter
Enums to identify who created this track and which properties does it have.
@ NumberOfTrackFitters
maximum number of enums

◆ trackLink()

const ElementLink< TrackCollection > & xAOD::TrackParticle_v1::trackLink ( ) const

Returns a link (which can be invalid) to the Trk::Track which was used to make this TrackParticle.

The function will return an invalid ElementLink in case nothing was set for it yet.

Note
This is only available in Athena.

This is to avoid users having to always check both for the decoration being available, and the link being valid.

Returns
An element link to the parent Trk::Track of this track particle

Definition at line 730 of file TrackParticle_v1.cxx.

730 {
731
732 // The accessor:
733 static const ConstAccessor< ElementLink< TrackCollection > > acc( "trackLink" );
734
735 // Check if one of them is available:
736 if( acc.isAvailable( *this ) ) {
737 return acc( *this );
738 }
739
740 // If no Trk::Track link was not set (yet), return a dummy object:
741 static const ElementLink< TrackCollection > dummy;
742 return dummy;
743 }

◆ trackParameterCovarianceMatrix()

xAOD::ParametersCovMatrix_t xAOD::TrackParticle_v1::trackParameterCovarianceMatrix ( unsigned int index) const

Returns the TrackParticleCovMatrix_t (covariance matrix) at 'index', which corresponds to the parameters at the same index.

Definition at line 579 of file TrackParticle_v1.cxx.

580 {
581 static const Accessor< std::vector<float> > acc( "trackParameterCovarianceMatrices" );
582 unsigned int offset = index*15;
583 // copy the correct values into the temp matrix
584 xAOD::ParametersCovMatrix_t tmp;
585 std::vector<float>::const_iterator it = acc(*this).begin()+offset;
586 Amg::expand(it,it+15,tmp);
587 return tmp;
588 }

◆ trackParameters()

const CurvilinearParameters_t xAOD::TrackParticle_v1::trackParameters ( unsigned int index) const

Returns the track parameter vector at 'index'.

Definition at line 512 of file TrackParticle_v1.cxx.

512 {
513 CurvilinearParameters_t tmp;
514 tmp << parameterX(index),parameterY(index),parameterZ(index),
515 parameterPX(index),parameterPY(index),parameterPZ(index);
516 return tmp;
517 }

◆ trackProperties()

TrackProperties xAOD::TrackParticle_v1::trackProperties ( ) const

Access methods for track properties, which returns 'true' if a logical AND of the parameter 'proprty' and the stored properties returns true.

i.e. you do:

TrackProperties testProperty;
testProperty.set(SOMEPROPERTY);
if (trackParticle.trackProperties(testProperty)) doSomething();
Todo
  • fix the above (or make something nicer)

◆ type()

Type::ObjectType xAOD::TrackParticle_v1::type ( ) const
finaloverridevirtual

The type of the object as a simple enumeration.

Implements xAOD::IParticle.

Definition at line 139 of file TrackParticle_v1.cxx.

139 {
140 return Type::TrackParticle;
141 }

◆ vx()

float xAOD::TrackParticle_v1::vx ( ) const

The x origin for the parameters.

◆ vy()

float xAOD::TrackParticle_v1::vy ( ) const

The y origin for the parameters.

◆ vz()

float xAOD::TrackParticle_v1::vz ( ) const

The z origin for the parameters.

◆ z0()

float xAOD::TrackParticle_v1::z0 ( ) const

Returns the \(z_0\) parameter.

Member Data Documentation

◆ COVMATRIX_OFFDIAG_VEC_COMPR_SIZE

const std::size_t xAOD::TrackParticle_v1::COVMATRIX_OFFDIAG_VEC_COMPR_SIZE = TrackingDetails::COVMATRIX_OFFDIAG_VEC_COMPR_SIZE
staticprivate

Definition at line 350 of file TrackParticle_v1.h.

◆ m_perigeeParameters

CxxUtils::CachedValue<Trk::Perigee> xAOD::TrackParticle_v1::m_perigeeParameters
private

Cached MeasuredPerigee, built from this object.

Note
This is only available in Athena.

Definition at line 369 of file TrackParticle_v1.h.


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