ATLAS Offline Software
Loading...
Searching...
No Matches
TrackSummary_v1.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4#ifndef XAODTRACKING_VERSIONS_TRACKSUMMARY_V1_H
5#define XAODTRACKING_VERSIONS_TRACKSUMMARY_V1_H
6#include <cstdint>
7#include <vector>
8#include "AthLinks/ElementLink.h"
12
13
14namespace xAOD
15{
19
21 {
22 private:
24 public:
25 TrackSummary_v1() = default;
29
30 template <std::size_t measdim = 6>
31 Eigen::Map<const Eigen::Matrix<double, measdim, 1>> paramsEigen() const
32 {
33 return Eigen::Map<const Eigen::Matrix<double, measdim, 1>>{s_paramsAcc(*this).data()};
34 }
35
38 template <std::size_t measdim = 6>
39 Eigen::Map<Eigen::Matrix<double, measdim, 1>> paramsEigen()
40 {
41 return Eigen::Map<Eigen::Matrix<double, measdim, 1>>{s_paramsAcc(*this).data()};
42 }
47 const std::vector<double> &params() const;
51 void setParams(const std::vector<double> &m);
52
56 template <std::size_t measdim = 6>
57 Eigen::Map<const Eigen::Matrix<double, measdim, measdim>> covParamsEigen() const
58 {
59 return Eigen::Map<const Eigen::Matrix<double, measdim, measdim>>{s_covParamsAcc(*this).data()};
60 }
61
65 template <std::size_t measdim = 6>
66 Eigen::Map<Eigen::Matrix<double, measdim, measdim>> covParamsEigen()
67 {
68 return Eigen::Map<Eigen::Matrix<double, measdim, measdim>>{s_covParamsAcc(*this).data()};
69 }
70
74 const std::vector<double> &covParams() const;
78 void setCovParams(const std::vector<double> &m);
79
83 unsigned int nMeasurements() const;
87 void setnMeasurements(unsigned int m);
91 const unsigned int* nMeasurementsPtr() const;
92 unsigned int* nMeasurementsPtr();
93
97 unsigned int nHoles() const;
101 void setnHoles(unsigned int m);
105 const unsigned int* nHolesPtr() const;
106 unsigned int* nHolesPtr();
107
108
109
113 float chi2f() const;
117 void setChi2f(float m);
121 const float* chi2fPtr() const;
122 float* chi2fPtr();
123
127 unsigned int ndf() const;
131 void setNdf(unsigned int m);
135 const unsigned int* ndfPtr() const;
136 unsigned int* ndfPtr();
137
138
142 unsigned int nOutliers() const;
146 void setnOutliers(unsigned int m);
150 const unsigned int* nOutliersPtr() const;
151 unsigned int* nOutliersPtr();
152
156 unsigned int nSharedHits() const;
160 void setnSharedHits(unsigned int m);
164 const unsigned int* nSharedHitsPtr() const;
165 unsigned int* nSharedHitsPtr();
166
167
171 unsigned int tipIndex() const;
172
177 void setTipIndex(unsigned int);
178
182 const unsigned int* tipIndexPtr() const;
183 unsigned int* tipIndexPtr();
184
189 unsigned int stemIndex() const;
190
195 void setStemIndex(unsigned int);
196
200 const unsigned int* stemIndexPtr() const;
201 unsigned int* stemIndexPtr();
202
203
207 unsigned int surfaceIndex() const;
208
212 void setSurfaceIndex(unsigned int);
213
214
215
221
225 void resize(size_t sz = 6);
226
230 size_t size() const;
231 };
232}
233#endif
Base class for elements of a container that can have aux data.
static Double_t sz
Base class for elements of a container that can have aux data.
Definition AuxElement.h:483
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
void resize(size_t sz=6)
resize internal arrays to store params (to capacity sz) & convariances (to capacity sz x sz)
unsigned int nSharedHits() const
access nSharedHits
void setnSharedHits(unsigned int m)
access set nSharedHits
const unsigned int * nHolesPtr() const
pointers API needed by MTJ
const unsigned int * tipIndexPtr() const
pointers API needed by MTJ
unsigned int nMeasurements() const
access nMeasurements
void setStemIndex(unsigned int)
Set the stem index.
const unsigned int * nMeasurementsPtr() const
pointers API needed by MTJ
float chi2f() const
access chi2
unsigned int * nMeasurementsPtr()
const unsigned int * stemIndexPtr() const
pointers API needed by MTJ
Eigen::Map< const Eigen::Matrix< double, measdim, 1 > > paramsEigen() const
access track backend vector of const element
static const SG::AuxElement::Accessor< std::vector< double > > s_covParamsAcc
void setChi2f(float m)
access set chi2
void setnMeasurements(unsigned int m)
access set nMeasurements
Eigen::Map< Eigen::Matrix< double, measdim, 1 > > paramsEigen()
access parameters of non const element
const unsigned int * nSharedHitsPtr() const
pointers API needed by MTJ
unsigned int * nHolesPtr()
size_t size() const
retrieve the size of the internal vectors for the data summary
Eigen::Map< const Eigen::Matrix< double, measdim, measdim > > covParamsEigen() const
access track covariance matrix (flattened, rows layout) of const element
const std::vector< double > & covParams() const
access track covariance as plain vector
unsigned int tipIndex() const
index of the tip of the TrackStates linked list in MultiTrajectory
const unsigned int * nOutliersPtr() const
pointers API needed by MTJ
const unsigned int * ndfPtr() const
pointers API needed by MTJ
void setTipIndex(unsigned int)
Set the tip index.
void setNdf(unsigned int m)
access set ndf
unsigned int * tipIndexPtr()
void setCovParams(const std::vector< double > &m)
access set covariance from plain vector
void setSurfaceIndex(unsigned int)
Set the index in surface container.
unsigned int surfaceIndex() const
index of the surfaces in the surfaces collection
unsigned int stemIndex() const
index of the stem of the TrackStates linked list in MultiTrajectory steemIndex is needed reverse look...
unsigned int nHoles() const
access nHoles
const std::vector< double > & params() const
access track parameters as plain vector
unsigned int nOutliers() const
access nOutliers
unsigned int * nOutliersPtr()
unsigned int * ndfPtr()
Eigen::Map< Eigen::Matrix< double, measdim, measdim > > covParamsEigen()
access track covariance matrix (flattened, rows layout)
unsigned int ndf() const
access ndf
static const SG::AuxElement::Accessor< std::vector< double > > s_paramsAcc
unsigned int * nSharedHitsPtr()
void setParticleHypothesis(const uint8_t &)
const float * chi2fPtr() const
pointers API needed by MTJ
void setnOutliers(unsigned int m)
access set nOutliers
const uint8_t & particleHypothesis() const
particle hypothesis access
void setParams(const std::vector< double > &m)
access set parameters from plain vector
unsigned int * stemIndexPtr()
void setnHoles(unsigned int m)
access set nHoles
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.