ATLAS Offline Software
FitMatrices.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // class FitMatrices
7 // Storage and manipulation of matrices during track fitting
8 // (note the actual matrices are structs (FitMatrix.h) to give faster
9 // execution)
10 //
11 // Given matrix of measurement derivatives wrt fit parameters (DerivativeMatrix
12 // DM) and vector of differences between measurements and fitted trajectory,
13 // solve for:
14 // parameter weight matrix = (covariance)-1 = DMtranspose.DM
15 // parameter change = (DMtranspose.DM)-1 * (DMtranspose.differences)
16 //
17 // NOTE:
18 // covariances, derivatives etc, use d0, z0, phi, cot(theta), qOverPt as first
19 // 5 parameters distinguish: full covariance with all parameters:
20 // includes misalignments, scattering and energy loss
21 // 5*5 final covariance with possible external contribution
22 // representing leading material and field gradient effects
23 //
25 
26 #ifndef TRKIPATFITTERUTILS_FITMATRICES_H
27 #define TRKIPATFITTERUTILS_FITMATRICES_H
28 
29 #include <cmath>
30 #include <list>
31 
33 #include "GaudiKernel/MsgStream.h"
36 
37 namespace Trk {
38 class FitMeasurement;
39 class FitParameters;
40 
41 class FitMatrices {
42  public:
43  FitMatrices(bool constrainedAlignmentEffects);
44 
45  // copy, assignment deleted
46  FitMatrices(const FitMatrices&) = delete;
47  FitMatrices& operator=(const FitMatrices&) = delete;
48  //move defaulted
49  FitMatrices(FitMatrices&&) = default;
51  //dtor default
53 
54  // debugging aid: check 'smart' pointers
55  void checkPointers(MsgStream& log) const;
56 
57  // change to chiSquared
58  static double chiSquaredChange(void);
59 
60  // accessor to DerivativeMatrix (eigen)
62 
63  // accessor to final covariance (5*5)
64  // includes leading material effects on construction, but
65  // field gradient effects have to be added by MeasurementProcessor as
66  // 'extrapolation aware'
68 
69  // produce full covariance (including scatterer parameters) i.e. invert weight
70  // matrix NOTE: this does not contain the external errors which will be
71  // included in the finalCovariance
72  // or the leading material contribution
73  const Amg::MatrixX* fullCovariance(void);
74 
75  // include external errors
76  // void externalMomentumUncertainty (double
77  // fractionalErrorSq);
78 
79  // number of degrees of freedom
80  int numberDoF(void) const;
81 
82  // number of driftCircles to fit
83  int numberDriftCircles(void) const;
84 
85  // chiSquared contribution from perigee measurement
86  double perigeeChiSquared(void);
87 
88  // debugging aids (n.b. using std::cout)
89  void printDerivativeMatrix(void);
90  void printWeightMatrix(void);
91 
92  // remove leading+trailing zeroes from smart pointers
93  void refinePointers(void);
94 
95  // release memory (eigen)
96  void releaseMemory(void);
97 
98  // initialize matrices - set appropriate dimensions for a given set of
99  // measurements
100  int setDimensions(std::vector<FitMeasurement*>& measurements,
102 
103  // solve matrix equations for parameters change
104  bool solveEquations(void);
105 
106  // tell fit to use measured perigee parameters
107  void usePerigee(const FitMeasurement& measurement);
108 
109  private:
110  // add perigee measurement
111  void addPerigeeMeasurement(void);
112  // fix for momentum singularity
113  void avoidMomentumSingularity(void); // using Eigen
114 
121  std::vector<int> m_firstRowForParameter;
123  std::vector<int> m_lastRowForParameter;
125  std::vector<FitMeasurement*>* m_measurements; // not owning ptr
130  const Amg::VectorX* m_perigee; // not owning ptr
132  const Amg::MatrixX* m_perigeeWeight; // not owning ptr
133  std::vector<double> m_residuals;
134  int m_rowsDM;
138 };
139 
140 //<<<<<< INLINE MEMBER FUNCTIONS >>>>>>
141 
143  return m_derivativeMatrix;
144 }
145 
147  return m_finalCovariance;
148 }
149 
150 inline int FitMatrices::numberDoF(void) const {
151  return m_numberDoF;
152 }
153 
154 inline int FitMatrices::numberDriftCircles(void) const {
155  return m_numberDriftCircles;
156 }
157 
158 } // namespace Trk
159 
160 #endif // TRKIPATFITTERUTILS_FITMATRICES_H
Trk::fitMatrix
Definition: FitMatrix.h:13
Trk::FitMatrices::m_weightedDifference
Amg::VectorX m_weightedDifference
Definition: FitMatrices.h:137
Amg::VectorX
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
Definition: EventPrimitives.h:32
Trk::FitMatrices::finalCovariance
Amg::MatrixX & finalCovariance()
Definition: FitMatrices.h:146
Amg::MatrixX
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Definition: EventPrimitives.h:29
Trk::FitMatrices::m_perigeeDifference
Amg::MatrixX m_perigeeDifference
Definition: FitMatrices.h:131
Trk::FitMatrices::avoidMomentumSingularity
void avoidMomentumSingularity(void)
Definition: FitMatrices.cxx:634
Trk::FitMatrices::m_weight
Amg::MatrixX m_weight
Definition: FitMatrices.h:136
Trk::FitMatrices::perigeeChiSquared
double perigeeChiSquared(void)
Definition: FitMatrices.cxx:162
Trk::FitMatrices::usePerigee
void usePerigee(const FitMeasurement &measurement)
Definition: FitMatrices.cxx:610
Trk::FitMatrices::numberDoF
int numberDoF(void) const
Definition: FitMatrices.h:150
Trk::FitMatrices::m_derivativeMatrix
Amg::MatrixX m_derivativeMatrix
Definition: FitMatrices.h:119
Trk::FitMatrices::m_rowsDM
int m_rowsDM
Definition: FitMatrices.h:134
Trk::FitMatrices::~FitMatrices
~FitMatrices()=default
Trk::FitMatrices::m_perigee
const Amg::VectorX * m_perigee
Definition: FitMatrices.h:130
Trk::FitMatrices::m_fitMatrix
fitMatrix m_fitMatrix
Definition: FitMatrices.h:115
Trk::FitMatrices::numberDriftCircles
int numberDriftCircles(void) const
Definition: FitMatrices.h:154
Trk::FitMatrices::refinePointers
void refinePointers(void)
Definition: FitMatrices.cxx:273
FitMatrix.h
Trk::FitMatrices::operator=
FitMatrices & operator=(const FitMatrices &)=delete
Trk::FitMatrices::m_numberDriftCircles
int m_numberDriftCircles
Definition: FitMatrices.h:127
Trk::FitMatrices::printWeightMatrix
void printWeightMatrix(void)
Definition: FitMatrices.cxx:252
GeoPrimitives.h
Trk::FitMatrices::chiSquaredChange
static double chiSquaredChange(void)
Definition: FitMatrices.cxx:81
Trk::FitMatrices::m_lastRowForParameter
std::vector< int > m_lastRowForParameter
Definition: FitMatrices.h:123
Trk::FitMeasurement
Definition: FitMeasurement.h:40
Trk::FitMatrices::FitMatrices
FitMatrices(bool constrainedAlignmentEffects)
Definition: FitMatrices.cxx:41
Trk::FitMatrices::m_usePerigee
bool m_usePerigee
Definition: FitMatrices.h:135
Trk::FitMatrices::printDerivativeMatrix
void printDerivativeMatrix(void)
Definition: FitMatrices.cxx:168
Trk::FitMatrices
Definition: FitMatrices.h:41
Trk::FitMatrices::m_finalCovariance
Amg::MatrixX m_finalCovariance
Definition: FitMatrices.h:120
Trk::FitMatrices::m_matrixFromCLHEP
bool m_matrixFromCLHEP
Definition: FitMatrices.h:124
EventPrimitives.h
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::FitMatrices::m_numberPerigee
int m_numberPerigee
Definition: FitMatrices.h:128
Trk::FitMatrices::m_residuals
std::vector< double > m_residuals
Definition: FitMatrices.h:133
Trk::FitMatrices::FitMatrices
FitMatrices(FitMatrices &&)=default
Trk::FitMatrices::FitMatrices
FitMatrices(const FitMatrices &)=delete
Trk::FitMatrices::m_largePhiWeight
double m_largePhiWeight
Definition: FitMatrices.h:122
Trk::FitMatrices::m_covariance
Amg::MatrixX m_covariance
Definition: FitMatrices.h:118
Trk::FitMatrices::derivativeMatrix
Amg::MatrixX & derivativeMatrix()
Definition: FitMatrices.h:142
Trk::FitMatrices::m_parameters
FitParameters * m_parameters
Definition: FitMatrices.h:129
Trk::FitMatrices::checkPointers
void checkPointers(MsgStream &log) const
Definition: FitMatrices.cxx:58
Trk::FitMatrices::addPerigeeMeasurement
void addPerigeeMeasurement(void)
Definition: FitMatrices.cxx:619
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
Trk::FitMatrices::fullCovariance
const Amg::MatrixX * fullCovariance(void)
Definition: FitMatrices.cxx:88
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
Trk::FitMatrices::m_columnsDM
int m_columnsDM
Definition: FitMatrices.h:116
Trk::FitParameters
Definition: FitParameters.h:29
python.CaloScaleNoiseConfig.default
default
Definition: CaloScaleNoiseConfig.py:79
Trk::FitMatrices::m_constrainedAlignmentEffects
bool m_constrainedAlignmentEffects
Definition: FitMatrices.h:117
Trk::FitMatrices::m_measurements
std::vector< FitMeasurement * > * m_measurements
Definition: FitMatrices.h:125
Trk::FitMatrices::m_numberDoF
int m_numberDoF
Definition: FitMatrices.h:126
Trk::FitMatrices::operator=
FitMatrices & operator=(FitMatrices &&)=default
Trk::FitMatrices::setDimensions
int setDimensions(std::vector< FitMeasurement * > &measurements, FitParameters *parameters)
Definition: FitMatrices.cxx:297
Trk::FitMatrices::releaseMemory
void releaseMemory(void)
Definition: FitMatrices.cxx:291
Trk::FitMatrices::m_firstRowForParameter
std::vector< int > m_firstRowForParameter
Definition: FitMatrices.h:121
Trk::FitMatrices::solveEquations
bool solveEquations(void)
Definition: FitMatrices.cxx:575
Trk::FitMatrices::m_perigeeWeight
const Amg::MatrixX * m_perigeeWeight
Definition: FitMatrices.h:132