ATLAS Offline Software
KalmanUpdatorSMatrix.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // KalmanUpdatorSMatrix.h
7 // Header file for Kalman Fitter Updator
9 // (c) ATLAS Detector software
11 // Markus.Elsing@cern.ch, Wolfgang.Liebig@cern.ch
13 
14 #ifndef TRK_KALMANUPDATORSMATRIX_H
15 #define TRK_KALMANUPDATORSMATRIX_H
16 
17 
22 // Eigen/Amg
24 // ROOT SMatrix lin alg: for internal calculations
25 #include "Math/SMatrix.h"
26 #include "Math/SVector.h"
27 #include <cmath>
28 #include <string_view>
29 
30 typedef ROOT::Math::SMatrix<double, 1, 1, ROOT::Math::MatRepSym<double, 1>>
32 typedef ROOT::Math::SMatrix<double, 2, 2, ROOT::Math::MatRepSym<double, 2>>
34 typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3>>
36 typedef ROOT::Math::SMatrix<double, 4, 4, ROOT::Math::MatRepSym<double, 4>>
38 typedef ROOT::Math::SMatrix<double, 5, 5, ROOT::Math::MatRepSym<double, 5>>
40 typedef ROOT::Math::SMatrix<double, 2, 2, ROOT::Math::MatRepStd<double, 2, 2>>
42 typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepStd<double, 3, 3>>
44 typedef ROOT::Math::SMatrix<double, 4, 4, ROOT::Math::MatRepStd<double, 4, 4>>
46 typedef ROOT::Math::SMatrix<double, 5, 5, ROOT::Math::MatRepStd<double, 5, 5>>
48 typedef ROOT::Math::SVector<double, 2> SParVector2;
49 typedef ROOT::Math::SVector<double, 3> SParVector3;
50 typedef ROOT::Math::SVector<double, 4> SParVector4;
51 typedef ROOT::Math::SVector<double, 5> SParVector5;
52 
53 namespace Trk {
54 
56 {
59 };
60 
75  : virtual public IUpdator
76  , public AthAlgTool
77 {
78 public:
80  KalmanUpdatorSMatrix(const std::string&,
81  const std::string&,
82  const IInterface*);
84 
86  virtual StatusCode initialize() override;
88  virtual StatusCode finalize() override;
89 
92  // fails: @copydoc Trk::IUpdator::addToState(const TrackParameters&, const
93  // Amg::Vector2D&, const Amg::MatrixX&)
94  virtual std::unique_ptr<TrackParameters> addToState(
95  const TrackParameters&,
96  const Amg::Vector2D&,
97  const Amg::MatrixX&) const override final;
100  virtual std::unique_ptr<TrackParameters> addToState(
101  const TrackParameters&,
102  const LocalParameters&,
103  const Amg::MatrixX&) const override final;
106  virtual std::unique_ptr<TrackParameters> addToState(
107  const TrackParameters&,
108  const Amg::Vector2D&,
109  const Amg::MatrixX&,
110  FitQualityOnSurface*&) const override final;
113  virtual std::unique_ptr<TrackParameters> addToState(
114  const TrackParameters&,
115  const LocalParameters&,
116  const Amg::MatrixX&,
117  FitQualityOnSurface*&) const override final;
118 
121  virtual std::unique_ptr<TrackParameters> removeFromState(
122  const TrackParameters&,
123  const Amg::Vector2D&,
124  const Amg::MatrixX&) const override final;
127  virtual std::unique_ptr<TrackParameters> removeFromState(
128  const TrackParameters&,
129  const LocalParameters&,
130  const Amg::MatrixX&) const override final;
133  virtual std::unique_ptr<TrackParameters> removeFromState(
134  const TrackParameters&,
135  const Amg::Vector2D&,
136  const Amg::MatrixX&,
137  FitQualityOnSurface*&) const override final;
140  virtual std::unique_ptr<TrackParameters> removeFromState(
141  const TrackParameters&,
142  const LocalParameters&,
143  const Amg::MatrixX&,
144  FitQualityOnSurface*&) const override final;
145 
153  virtual std::unique_ptr<TrackParameters> combineStates(
154  const TrackParameters&,
155  const TrackParameters&) const override final;
162  virtual std::unique_ptr<TrackParameters> combineStates(
163  const TrackParameters&,
164  const TrackParameters&,
165  FitQualityOnSurface*&) const override final;
166 
170  const TrackParameters&,
171  const Amg::Vector2D&,
172  const Amg::MatrixX&) const override final;
176  const TrackParameters&,
177  const LocalParameters&,
178  const Amg::MatrixX&) const override final;
182  const TrackParameters&,
183  const Amg::Vector2D&,
184  const Amg::MatrixX&) const override final;
188  const TrackParameters&,
189  const LocalParameters&,
190  const Amg::MatrixX&) const override final;
195  const TrackParameters&,
196  const TrackParameters&) const override final;
199  const AmgVector(5) &,
200  const AmgSymMatrix(5) &,
201  const Amg::VectorX&,
202  const Amg::MatrixX&,
203  int,
205  bool) const override final
206  {
207  return nullptr;
208  }
209 
211  virtual std::vector<double> initialErrors() const override final;
212 
213 private:
218  const Amg::MatrixX&,
219  const int,
221  bool) const;
224  std::unique_ptr<TrackParameters> calculateFilterStep_1D(
228  double,
229  int,
230  const Amg::MatrixX&,
231  const int,
233  bool) const;
236  std::unique_ptr<TrackParameters> calculateFilterStep_2D(
241  int,
242  const Amg::MatrixX&,
243  const int,
245  bool) const;
248  std::unique_ptr<TrackParameters> calculateFilterStep_3D(
253  const Amg::MatrixX&,
254  const int,
256  bool) const;
259  std::unique_ptr<TrackParameters> calculateFilterStep_4D(
264  const Amg::MatrixX&,
265  const int,
267  bool) const;
271 
273  std::unique_ptr<TrackParameters> calculateFilterStep_5D(
278  const Amg::MatrixX&,
279  const int,
281  bool) const;
282 
286  std::unique_ptr<TrackParameters> convertToClonedTrackPars(
290  int,
291  bool,
292  std::string_view) const;
302  const AmgSymMatrix(5)&,
303  double,
304  double,
305  int,
306  int) const;
308  const AmgSymMatrix(5)&,
311  int,
312  int) const;
314  const AmgSymMatrix(5)&,
316  const AmgSymMatrix(5)&,
317  int) const;
319  const AmgSymMatrix(5)&,
320  const Amg::MatrixX&,
321  const Amg::MatrixX&,
322  int) const;
323 
329  static SCovMatrix2 projection_2D(const Amg::MatrixX&, int) ;
336 
337  // === note: any of the following log... method is only called if
338  // the msgstream level has been set appropriately.
340  void logStart(const std::string&, const TrackParameters&) const;
343  const Amg::VectorX&,
344  const Amg::MatrixX&) const;
346  void logGainForm(int,
351  void logResult(const std::string&,
352  const AmgVector(5) &,
353  const AmgSymMatrix(5) &) const;
354 
357 
359 
367 
370 
373 
376  SCovMatrix5&,
378 
379  std::vector<double> m_cov_stdvec;
384 
387 };
388 inline bool
390  const RangeCheckDef rcd) const
391 {
392  static const SParVector2 thetaMin(0.0, -M_PI);
393  return ((std::abs(V(Trk::phi)) <= M_PI) &&
394  (V(Trk::theta) >= thetaMin((int)rcd)) && (V(Trk::theta) <= M_PI));
395 }
396 
397 inline bool
399 {
400  return (V(Trk::theta) >= 0.0 && (V(Trk::theta) <= M_PI));
401 }
402 
403 } // end of namespace
404 
405 #endif // TRK_KALMANUPDATOR_H
Trk::KalmanUpdatorSMatrix::projection_2D
static SCovMatrix2 projection_2D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 2D rows&columns out of the full cov matrix.
Definition: KalmanUpdatorSMatrix.cxx:1193
Trk::LocalParameters
Definition: LocalParameters.h:98
SParVector3
ROOT::Math::SVector< double, 3 > SParVector3
Definition: KalmanUpdatorSMatrix.h:49
Trk::KalmanUpdatorSMatrix::fullStateFitQuality
virtual FitQualityOnSurface fullStateFitQuality(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
estimator for FitQuality on Surface from a full track state, that is a state which contains the curre...
Definition: KalmanUpdatorSMatrix.cxx:335
Trk::KalmanUpdatorSMatrix::~KalmanUpdatorSMatrix
virtual ~KalmanUpdatorSMatrix()
Trk::KalmanUpdatorSMatrix::updateParameterDifference
virtual std::pair< AmgVector(5), AmgSymMatrix(5)> * updateParameterDifference(const AmgVector(5) &, const AmgSymMatrix(5) &, const Amg::VectorX &, const Amg::MatrixX &, int, Trk::FitQualityOnSurface *&, bool) const override final
interface for reference-track KF, not implemented
Definition: KalmanUpdatorSMatrix.h:198
Amg::VectorX
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
Definition: EventPrimitives.h:30
Amg::MatrixX
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Definition: EventPrimitives.h:27
TrackParameters.h
Trk::KalmanUpdatorSMatrix::combineStates
virtual std::unique_ptr< TrackParameters > combineStates(const TrackParameters &, const TrackParameters &) const override final
trajectory state updator which combines two parts of a trajectory on a common surface.
Definition: KalmanUpdatorSMatrix.cxx:252
Trk::KalmanUpdatorSMatrix::makeChi2_2D
FitQualityOnSurface makeChi2_2D(const SParVector5 &, const AmgSymMatrix(5)&, const SParVector2 &, const SCovMatrix2 &, int, int) const
Definition: KalmanUpdatorSMatrix.cxx:1103
Trk::differentialCheck
@ differentialCheck
Definition: KalmanUpdatorSMatrix.h:58
Trk::KalmanUpdatorSMatrix::thetaWithinRange_5D
bool thetaWithinRange_5D(const SParVector5 &V) const
Test if theta angle is inside boundaries. No differential-check option.
Definition: KalmanUpdatorSMatrix.h:398
Trk::KalmanUpdatorSMatrix::projection_4D
static SCovMatrix4 projection_4D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 4D rows&columns out of the full cov matrix.
Definition: KalmanUpdatorSMatrix.cxx:1242
Amg::Vector2D
Eigen::Matrix< double, 2, 1 > Vector2D
Definition: GeoPrimitives.h:48
Trk::KalmanUpdatorSMatrix::projection_3D
static SCovMatrix3 projection_3D(const SCovMatrix5 &, int)
Avoid multiplications with sparse H matrices by cutting 3D rows&columns out of the full cov matrix.
Definition: KalmanUpdatorSMatrix.cxx:1225
Trk::KalmanUpdatorSMatrix::s_enumAccessor
static const ParamDefsAccessor s_enumAccessor
Definition: KalmanUpdatorSMatrix.h:386
Trk::KalmanUpdatorSMatrix::convertToClonedTrackPars
std::unique_ptr< TrackParameters > convertToClonedTrackPars(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, int, bool, std::string_view) const
Helper method to convert internal results from SMatrix to Eigen. *‍/.
Definition: KalmanUpdatorSMatrix.cxx:1159
Trk::KalmanUpdatorSMatrix::initialize
virtual StatusCode initialize() override
AlgTool initialisation.
Definition: KalmanUpdatorSMatrix.cxx:50
M_PI
#define M_PI
Definition: ActiveFraction.h:11
SCovMatrix3
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > SCovMatrix3
Definition: KalmanUpdatorSMatrix.h:35
Trk::KalmanUpdatorSMatrix::initialErrors
virtual std::vector< double > initialErrors() const override final
give back how updator is configured for inital covariances
Definition: KalmanUpdatorSMatrix.cxx:541
SParVector2
ROOT::Math::SVector< double, 2 > SParVector2
Definition: KalmanUpdatorSMatrix.h:48
SGenMatrix4
ROOT::Math::SMatrix< double, 4, 4, ROOT::Math::MatRepStd< double, 4, 4 > > SGenMatrix4
Definition: KalmanUpdatorSMatrix.h:45
Trk::KalmanUpdatorSMatrix::logInputCov
void logInputCov(const SCovMatrix5 &, const Amg::VectorX &, const Amg::MatrixX &) const
internal structuring: common logfile output of the inputs
Definition: KalmanUpdatorSMatrix.cxx:1375
Trk::KalmanUpdatorSMatrix::m_unitMatrix
SCovMatrix5 m_unitMatrix
avoid mem allocation at every call
Definition: KalmanUpdatorSMatrix.h:385
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:560
Trk::KalmanUpdatorSMatrix::logStart
void logStart(const std::string &, const TrackParameters &) const
internal structuring: debugging output for start of method.
Definition: KalmanUpdatorSMatrix.cxx:1365
ParamDefs.h
Trk::KalmanUpdatorSMatrix
Implementation of Trk::IUpdator based on gain formalism and SMatrix mathlib.
Definition: KalmanUpdatorSMatrix.h:77
IUpdator.h
Trk::KalmanUpdatorSMatrix::consistentParamDimensions
bool consistentParamDimensions(const LocalParameters &, int) const
method testing correct use of LocalParameters
Definition: KalmanUpdatorSMatrix.cxx:1310
SGenMatrix2
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepStd< double, 2, 2 > > SGenMatrix2
Definition: KalmanUpdatorSMatrix.h:41
Trk::AmgSymMatrix
AmgSymMatrix(5) &GXFTrackState
Definition: GXFTrackState.h:156
Trk::FitQualityOnSurface
Definition: FitQualityOnSurface.h:19
Trk::KalmanUpdatorSMatrix::calculateFilterStep_3D
std::unique_ptr< TrackParameters > calculateFilterStep_3D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
Definition: KalmanUpdatorSMatrix.cxx:812
SCovMatrix2
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > SCovMatrix2
Definition: KalmanUpdatorSMatrix.h:33
Trk::KalmanUpdatorSMatrix::m_cov_stdvec
std::vector< double > m_cov_stdvec
job options for initial cov values
Definition: KalmanUpdatorSMatrix.h:379
Trk::KalmanUpdatorSMatrix::getStartCov
bool getStartCov(SCovMatrix5 &, const TrackParameters &, const int) const
Helper method to transform Eigen cov matrix to SMatrix.
Definition: KalmanUpdatorSMatrix.cxx:1259
Trk::KalmanUpdatorSMatrix::m_cov0
SParVector5 m_cov0
initial cov values in SMatrix object
Definition: KalmanUpdatorSMatrix.h:380
Trk::theta
@ theta
Definition: ParamDefs.h:66
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
Trk::KalmanUpdatorSMatrix::correctThetaPhiRange_5D
bool correctThetaPhiRange_5D(SParVector5 &, SCovMatrix5 &, const RangeCheckDef) const
method correcting the calculated angles back to their defined ranges phi (-pi, pi) and theta (0,...
Definition: KalmanUpdatorSMatrix.cxx:1319
vector
Definition: MultiHisto.h:13
Trk::KalmanUpdatorSMatrix::calculateFilterStep_4D
std::unique_ptr< TrackParameters > calculateFilterStep_4D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
Definition: KalmanUpdatorSMatrix.cxx:904
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
AthAlgTool.h
SCovMatrix4
ROOT::Math::SMatrix< double, 4, 4, ROOT::Math::MatRepSym< double, 4 > > SCovMatrix4
Definition: KalmanUpdatorSMatrix.h:37
Trk::ParametersBase
Definition: ParametersBase.h:55
Trk::KalmanUpdatorSMatrix::addToState
virtual std::unique_ptr< TrackParameters > addToState(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
measurement updator for the KalmanFitter getting the meas't coord' from Amg::Vector2D (use eg with PR...
Definition: KalmanUpdatorSMatrix.cxx:76
Trk::KalmanUpdatorSMatrix::KalmanUpdatorSMatrix
KalmanUpdatorSMatrix(const std::string &, const std::string &, const IInterface *)
AlgTool standard constuctor.
Definition: KalmanUpdatorSMatrix.cxx:33
Trk::KalmanUpdatorSMatrix::logGainForm
void logGainForm(int, const SParVector5 &, const SCovMatrix5 &, const SGenMatrix5 &) const
internal structuring: common logfile output during calculation
Definition: KalmanUpdatorSMatrix.cxx:1401
SCovMatrix1
ROOT::Math::SMatrix< double, 1, 1, ROOT::Math::MatRepSym< double, 1 > > SCovMatrix1
Definition: KalmanUpdatorSMatrix.h:31
Trk::RangeCheckDef
RangeCheckDef
Definition: KalmanUpdatorSMatrix.h:56
Trk::KalmanUpdatorSMatrix::calculateFilterStep_5D
std::unique_ptr< TrackParameters > calculateFilterStep_5D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const SParVector5 &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
Definition: KalmanUpdatorSMatrix.cxx:997
EventPrimitives.h
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Amg
Definition of ATLAS Math & Geometry primitives (Amg)
Definition: AmgStringHelpers.h:19
Trk::KalmanUpdatorSMatrix::thetaPhiWithinRange_5D
bool thetaPhiWithinRange_5D(const SParVector5 &V, const RangeCheckDef rcd) const
Test if angles are inside boundaries.
Definition: KalmanUpdatorSMatrix.h:389
SCovMatrix5
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > SCovMatrix5
Definition: KalmanUpdatorSMatrix.h:39
private
#define private
Definition: DetDescrConditionsDict_dict_fixes.cxx:13
Trk::absoluteCheck
@ absoluteCheck
Definition: KalmanUpdatorSMatrix.h:57
Trk::IUpdator
Set of interfaces for methods operating on track states, mainly for Kalman filtering.
Definition: IUpdator.h:64
Trk::KalmanUpdatorSMatrix::logResult
void logResult(const std::string &, const AmgVector(5) &, const AmgSymMatrix(5) &) const
internal structuring: common logfile output after calculation
Definition: KalmanUpdatorSMatrix.cxx:1424
Trk::KalmanUpdatorSMatrix::makeChi2_5D
FitQualityOnSurface makeChi2_5D(const SParVector5 &, const AmgSymMatrix(5)&, const SParVector5 &, const AmgSymMatrix(5)&, int) const
Definition: KalmanUpdatorSMatrix.cxx:1130
SGenMatrix3
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepStd< double, 3, 3 > > SGenMatrix3
Definition: KalmanUpdatorSMatrix.h:43
Trk::KalmanUpdatorSMatrix::m_thetaGainDampingValue
float m_thetaGainDampingValue
Definition: KalmanUpdatorSMatrix.h:383
Trk::KalmanUpdatorSMatrix::calculateFilterStep_1D
std::unique_ptr< TrackParameters > calculateFilterStep_1D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, double, int, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
Definition: KalmanUpdatorSMatrix.cxx:590
Trk::ParamDefsAccessor
Definition: ParamDefs.h:92
Trk::KalmanUpdatorSMatrix::makeChi2_1D
FitQualityOnSurface makeChi2_1D(const SParVector5 &, const AmgSymMatrix(5)&, double, double, int, int) const
also the chi2 calculation and FitQuality object creation is combined in an extra method.
Definition: KalmanUpdatorSMatrix.cxx:1081
Trk::KalmanUpdatorSMatrix::predictedStateFitQuality
virtual FitQualityOnSurface predictedStateFitQuality(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
estimator for FitQuality on Surface from a predicted track state, that is a state which contains the ...
Definition: KalmanUpdatorSMatrix.cxx:427
Trk::KalmanUpdatorSMatrix::prepareFilterStep
std::unique_ptr< TrackParameters > prepareFilterStep(const TrackParameters &, const LocalParameters &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common code analysing the measurement's rank and calling the appropriate implementation for this rank...
Definition: KalmanUpdatorSMatrix.cxx:548
Trk::phi
@ phi
Definition: ParamDefs.h:75
Trk::KalmanUpdatorSMatrix::makeChi2Object
FitQualityOnSurface makeChi2Object(const Amg::VectorX &, const AmgSymMatrix(5)&, const Amg::MatrixX &, const Amg::MatrixX &, int) const
Definition: KalmanUpdatorSMatrix.cxx:1286
SParVector4
ROOT::Math::SVector< double, 4 > SParVector4
Definition: KalmanUpdatorSMatrix.h:50
AthAlgTool
Definition: AthAlgTool.h:26
Trk::KalmanUpdatorSMatrix::finalize
virtual StatusCode finalize() override
AlgTool termination.
Definition: KalmanUpdatorSMatrix.cxx:70
SGenMatrix5
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > SGenMatrix5
Definition: KalmanUpdatorSMatrix.h:47
Trk::KalmanUpdatorSMatrix::calculateFilterStep_2D
std::unique_ptr< TrackParameters > calculateFilterStep_2D(const TrackParameters &, const SParVector5 &, const SCovMatrix5 &, const SParVector2 &, int, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
common maths calculation code for all addToState and removeFromState versions which happen to be call...
Definition: KalmanUpdatorSMatrix.cxx:702
SParVector5
ROOT::Math::SVector< double, 5 > SParVector5
Definition: KalmanUpdatorSMatrix.h:51
Trk::KalmanUpdatorSMatrix::removeFromState
virtual std::unique_ptr< TrackParameters > removeFromState(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &) const override final
reverse update eg for track property analysis (unbiased residuals) getting the measurement coordinate...
Definition: KalmanUpdatorSMatrix.cxx:164
Trk::KalmanUpdatorSMatrix::m_useFruehwirth8a
bool m_useFruehwirth8a
job options controlling update formula for covariance matrix
Definition: KalmanUpdatorSMatrix.h:381