ATLAS Offline Software
Loading...
Searching...
No Matches
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
37namespace Trk {
38class FitMeasurement;
39class FitParameters;
40
42 public:
43 FitMatrices(bool constrainedAlignmentEffects);
44
45 // copy, assignment deleted
46 FitMatrices(const FitMatrices&) = delete;
48 //move defaulted
51 //dtor default
52 ~FitMatrices() = 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,
101 FitParameters* parameters);
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;
138};
139
140//<<<<<< INLINE MEMBER FUNCTIONS >>>>>>
141
145
149
150inline int FitMatrices::numberDoF(void) const {
151 return m_numberDoF;
152}
153
154inline int FitMatrices::numberDriftCircles(void) const {
156}
157
158} // namespace Trk
159
160#endif // TRKIPATFITTERUTILS_FITMATRICES_H
double perigeeChiSquared(void)
Amg::MatrixX m_perigeeDifference
void usePerigee(const FitMeasurement &measurement)
void refinePointers(void)
void releaseMemory(void)
Amg::MatrixX & finalCovariance()
int numberDriftCircles(void) const
std::vector< double > m_residuals
~FitMatrices()=default
FitMatrices & operator=(FitMatrices &&)=default
const Amg::MatrixX * fullCovariance(void)
void addPerigeeMeasurement(void)
Amg::MatrixX m_weight
void checkPointers(MsgStream &log) const
std::vector< FitMeasurement * > * m_measurements
const Amg::MatrixX * m_perigeeWeight
bool m_constrainedAlignmentEffects
Amg::VectorX m_weightedDifference
std::vector< int > m_firstRowForParameter
Amg::MatrixX m_covariance
FitMatrices & operator=(const FitMatrices &)=delete
int numberDoF(void) const
const Amg::VectorX * m_perigee
void printWeightMatrix(void)
std::vector< int > m_lastRowForParameter
FitMatrices(FitMatrices &&)=default
void avoidMomentumSingularity(void)
FitMatrices(bool constrainedAlignmentEffects)
Amg::MatrixX m_derivativeMatrix
Amg::MatrixX & derivativeMatrix()
FitMatrices(const FitMatrices &)=delete
bool solveEquations(void)
Amg::MatrixX m_finalCovariance
static double chiSquaredChange(void)
double m_largePhiWeight
FitParameters * m_parameters
int setDimensions(std::vector< FitMeasurement * > &measurements, FitParameters *parameters)
fitMatrix m_fitMatrix
void printDerivativeMatrix(void)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
Ensure that the ATLAS eigen extensions are properly loaded.