ATLAS Offline Software
KalmanUpdator.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 // KalmanUpdator.h
7 // Header file for Kalman Fitter Updator
9 // (c) ATLAS Detector software
11 // Markus.Elsing@cern.ch
13 
14 #ifndef TRK_KALMANUPDATOR_H
15 #define TRK_KALMANUPDATOR_H
16 
20 #include "GaudiKernel/MsgStream.h"
26 #include <cmath>
27 
28 namespace Trk {
29 
43 class KalmanUpdator final
44  : virtual public IUpdator
45  , public AthAlgTool
46 {
47 public:
49  KalmanUpdator(const std::string&, const std::string&, const IInterface*);
51 
53  virtual StatusCode initialize() override final;
55  virtual StatusCode finalize() override final;
56 
58  // from Amg::Vector2D (use eg with PRD)
59  virtual std::unique_ptr<TrackParameters> addToState(
61  const Amg::Vector2D&,
62  const Amg::MatrixX&) const override final;
64  // from LocalParameters (use eg with MeasurementBase, ROT)
65  virtual std::unique_ptr<TrackParameters> addToState(
68  const Amg::MatrixX&) const override final;
71  // of the state at the same time (Amg::Vector2D-version)
72  virtual std::unique_ptr<TrackParameters> addToState(
74  const Amg::Vector2D&,
75  const Amg::MatrixX&,
76  FitQualityOnSurface*&) const override final;
79  // of the state at the same time (LocalParameters-version)
80  virtual std::unique_ptr<TrackParameters> addToState(
83  const Amg::MatrixX&,
84  FitQualityOnSurface*&) const override final;
85 
87  // getting the measurement coordinates from the Amg::Vector2D class.
88  virtual std::unique_ptr<TrackParameters> removeFromState(
90  const Amg::Vector2D&,
91  const Amg::MatrixX&) const override final;
93  // getting the measurement coordinates from the LocalParameters class.
94  virtual std::unique_ptr<TrackParameters> removeFromState(
97  const Amg::MatrixX&) const override final;
99  // interface with Amg::Vector2D and FitQualityOnSurface.
100  virtual std::unique_ptr<TrackParameters> removeFromState(
102  const Amg::Vector2D&,
103  const Amg::MatrixX&,
104  FitQualityOnSurface*&) const override final;
106  // the interface with LocalParameters and FitQualityOnSurface.
107  virtual std::unique_ptr<TrackParameters> removeFromState(
110  const Amg::MatrixX&,
111  FitQualityOnSurface*&) const override final;
112 
119  virtual std::unique_ptr<TrackParameters> combineStates(
121  const TrackParameters&) const override final;
128  virtual std::unique_ptr<TrackParameters> combineStates(
131  FitQualityOnSurface*&) const override final;
132 
135  // which contains the current hit (expressed as Amg::Vector2D).
138  const Amg::Vector2D&,
139  const Amg::MatrixX&) const override final;
140 
142  // a state which contains the current hit (expressed as LocalParameters).
146  const Amg::MatrixX&) const override final;
149  // which contains the current hit (expressed as Amg::Vector2D).
152  const Amg::Vector2D&,
153  const Amg::MatrixX&) const override final;
156  // which contains the current hit (expressed as LocalParameters).
160  const Amg::MatrixX&) const override final;
163  // the parameters of another trajectory part extrapolated to the common
164  // surface.
167  const TrackParameters&) const override final;
170  const AmgVector(5) &,
171  const AmgSymMatrix(5) &,
172  const Amg::VectorX&,
173  const Amg::MatrixX&,
174  int,
176  bool) const override final
177  {
178  return nullptr;
179  }
180 
182  virtual std::vector<double> initialErrors() const override final;
183 
184 private:
188  const Amg::Vector2D&,
189  const Amg::MatrixX&,
190  const int,
192  bool) const;
197  const Amg::MatrixX&,
198  const int,
200  bool) const;
208  const Amg::MatrixX&,
209  const Amg::MatrixX&,
210  const int,
211  const int) const;
212 
215 
217  void logStart(const std::string&, const TrackParameters&) const;
219  void logInputCov(const Amg::MatrixX&,
220  const Amg::VectorX&,
221  const Amg::MatrixX&) const;
223  void logGainForm(int,
224  const Amg::VectorX&,
225  const Amg::MatrixX&,
226  const Amg::MatrixX&,
227  const Amg::MatrixX&) const;
229  void logResult(const std::string&,
230  const Amg::VectorX&,
231  const Amg::MatrixX&) const;
232 
236  bool thetaPhiWithinRange(const Amg::VectorX&, const int key = 31) const;
242  AmgSymMatrix(5) &,
243  const bool isDifference = false,
244  const int key = 31) const;
246  Amg::MatrixX&,
247  const bool isDifference = false,
248  const int key = 31) const;
249 
250  std::vector<double> m_cov0;
255 };
256 
257 inline bool
259 {
260  if (!(key & 4 || key & 8))
261  return true; // in case no angles measured.
262  if (key == 31)
263  return ((std::abs(V[Trk::phi]) <= M_PI) && (V[Trk::theta] >= 0.0) &&
264  (V[Trk::theta] <= M_PI));
265  else { // if vector is compressed (i.e. localParameters) need to extract
266  // phi,theta first.
267  bool okay = true;
268  if (key & 4) { // phi is being measured
269  int jphi = 0;
270  for (int itag = 0, ipos = 1; itag < Trk::phi; ++itag, ipos *= 2)
271  if (key & ipos)
272  ++jphi;
273  okay = okay && (std::abs(V[jphi]) <= M_PI);
274  }
275  if (key & 8) { // theta is being measured
276  int jtheta = 0;
277  for (int itag = 0, ipos = 1; itag <= Trk::theta; ++itag, ipos *= 2)
278  if (key & ipos)
279  ++jtheta;
280  okay = okay && (std::abs(V[jtheta] - M_PI_2) <= M_PI_2);
281  }
282  return okay;
283  }
284 }
285 
286 // for speed reasons make two separate inline functions
287 // phi differences should be smaller than pi (else go other way round) => same
288 // as absolute phi value. theta differences should be smaller than pi but can be
289 // negative => other constraint than absolute theta.
290 inline bool
292  const int key) const
293 {
294  if (!(key & 4 || key & 8))
295  return true; // in case no angles measured.
296  if (key == 31)
297  return ((std::abs(V[Trk::phi]) <= M_PI) && (V[Trk::theta] >= -M_PI) &&
298  (V[Trk::theta] <= M_PI));
299 
300  else { // if vector is compressed (i.e. localParameters) need to extract
301  // phi,theta first.
302  bool okay = true;
303  if (key & 4) { // phi is being measured
304  int jphi = 0;
305  for (int itag = 0, ipos = 1; itag < Trk::phi; ++itag, ipos *= 2)
306  if (key & ipos)
307  ++jphi;
308  okay = okay && (std::abs(V[jphi]) <= M_PI);
309  }
310  if (key & 8) { // theta is being measured
311  int jtheta = 0;
312  for (int itag = 0, ipos = 1; itag <= Trk::theta; ++itag, ipos *= 2)
313  if (key & ipos)
314  ++jtheta;
315  okay = okay && (std::abs(V[jtheta]) <= M_PI);
316  }
317  return okay;
318  }
319 }
320 
321 inline FitQualityOnSurface
323  const Amg::MatrixX& covTrk,
324  const Amg::MatrixX& covRio,
325  const int sign,
326  const int key) const
327 { // sign: -1 = updated, +1 = predicted parameters.
328  Amg::MatrixX R = covRio + sign * projection(covTrk, key); // .similarity(H);
329  if (R.determinant() == 0) {
330  ATH_MSG_DEBUG("matrix inversion failed");
331  return FitQualityOnSurface(0.0, (int)covRio.cols());
332  }
333  // get chi2 = r.T() * R^-1 * r
334  const double chiSquared = Amg::chi2(R.inverse(), residual);
335  ATH_MSG_VERBOSE("-U- fitQuality of " << (sign > 0 ? "predicted" : "updated")
336  << " state, chi2 :" << chiSquared
337  << " / ndof= " << covRio.cols());
338  // return the FitQualityOnSurface object
339  return FitQualityOnSurface(chiSquared, int(covRio.cols()));
340 }
341 
342 } // end of namespace
343 
344 #endif // TRK_KALMANUPDATOR_CLHEP_H
Trk::LocalParameters
Definition: LocalParameters.h:98
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
Trk::KalmanUpdator::diffThetaPhiWithinRange
bool diffThetaPhiWithinRange(const Amg::VectorX &, const int key=31) const
tests if ranges of phi (-pi, pi) and theta (0, pi) residuals are correct *‍/
Definition: KalmanUpdator.h:291
Trk::KalmanUpdator::initialize
virtual StatusCode initialize() override final
AlgTool initialisation.
Definition: KalmanUpdator.cxx:41
Trk::KalmanUpdator::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: KalmanUpdator.h:169
Trk::KalmanUpdator::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: KalmanUpdator.cxx:152
ClusterSeg::residual
@ residual
Definition: ClusterNtuple.h:20
Amg::Vector2D
Eigen::Matrix< double, 2, 1 > Vector2D
Definition: GeoPrimitives.h:48
Trk::ProjectionMatricesSet
the matrices to access the variably-dimensioned local parameters and map them to the defined five tra...
Definition: ProjectionMatricesSet.h:29
M_PI
#define M_PI
Definition: ActiveFraction.h:11
Trk::KalmanUpdator::finalize
virtual StatusCode finalize() override final
AlgTool termination.
Definition: KalmanUpdator.cxx:53
Trk::KalmanUpdator::~KalmanUpdator
~KalmanUpdator()
FitQualityOnSurface.h
Trk::KalmanUpdator::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
Definition: KalmanUpdator.cxx:337
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:560
ParamDefs.h
IUpdator.h
Trk::KalmanUpdator::makeChi2Object
FitQualityOnSurface makeChi2Object(const Amg::VectorX &, const Amg::MatrixX &, const Amg::MatrixX &, const int, const int) const
also the chi2 calculation and FitQuality object creation is combined in an extra method.
Definition: KalmanUpdator.h:322
Trk::AmgSymMatrix
AmgSymMatrix(5) &GXFTrackState
Definition: GXFTrackState.h:156
Trk::FitQualityOnSurface
Definition: FitQualityOnSurface.h:19
GeoPrimitives.h
Trk::KalmanUpdator::m_projectionMatrices
ProjectionMatricesSet m_projectionMatrices
get the correct projection matrix
Definition: KalmanUpdator.h:252
Trk::KalmanUpdator::KalmanUpdator
KalmanUpdator(const std::string &, const std::string &, const IInterface *)
AlgTool standard constuctor.
Definition: KalmanUpdator.cxx:24
Trk::theta
@ theta
Definition: ParamDefs.h:66
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
vector
Definition: MultiHisto.h:13
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
AthAlgTool.h
Trk::KalmanUpdator::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'
Definition: KalmanUpdator.cxx:59
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:107
Trk::KalmanUpdator::m_useFruehwirth8a
bool m_useFruehwirth8a
job option: formula for cov update
Definition: KalmanUpdator.h:253
Trk::KalmanUpdator::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
Definition: KalmanUpdator.cxx:273
Trk::ParametersBase
Definition: ParametersBase.h:55
Trk::KalmanUpdator::consistentParamDimensions
bool consistentParamDimensions(const LocalParameters &, int) const
method testing correct use of LocalParameters *‍/
Definition: KalmanUpdator.cxx:666
Trk::KalmanUpdator
Implementation of Trk::IUpdator based on gain formalism and Eigen.
Definition: KalmanUpdator.h:46
Trk::KalmanUpdator::logStart
void logStart(const std::string &, const TrackParameters &) const
internal structuring: debugging output for start of method.
Definition: KalmanUpdator.cxx:845
Trk::KalmanUpdator::logResult
void logResult(const std::string &, const Amg::VectorX &, const Amg::MatrixX &) const
internal structuring: common logfile output after calculation
Definition: KalmanUpdator.cxx:907
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::KalmanUpdator::thetaPhiWithinRange
bool thetaPhiWithinRange(const Amg::VectorX &, const int key=31) const
tests if ranges of abolute phi (-pi, pi) and theta (0, pi) are correct *‍/
Definition: KalmanUpdator.h:258
private
#define private
Definition: DetDescrConditionsDict_dict_fixes.cxx:13
Trk::IUpdator
Set of interfaces for methods operating on track states, mainly for Kalman filtering.
Definition: IUpdator.h:64
Trk::KalmanUpdator::logGainForm
void logGainForm(int, const Amg::VectorX &, const Amg::MatrixX &, const Amg::MatrixX &, const Amg::MatrixX &) const
internal structuring: common logfile output during calculation
Definition: KalmanUpdator.cxx:882
IDTPM::chiSquared
float chiSquared(const U &p)
Definition: TrackParametersHelper.h:126
EventPrimitivesCovarianceHelpers.h
Trk::KalmanUpdator::m_cov0
std::vector< double > m_cov0
job option: initial covariance matrix
Definition: KalmanUpdator.h:250
Trk::KalmanUpdator::correctThetaPhiRange
bool correctThetaPhiRange(Amg::VectorX &, AmgSymMatrix(5) &, const bool isDifference=false, const int key=31) const
brings phi/theta back into valid range using 2pi periodicity
Trk::KalmanUpdator::projection
Amg::MatrixX projection(const Amg::MatrixX &, const int) const
avoid CLHEP's empty math operations (H-matrix) by copying members out
Definition: KalmanUpdator.cxx:653
Trk::KalmanUpdator::initialErrors
virtual std::vector< double > initialErrors() const override final
gives back how updator is configured for inital covariances
Definition: KalmanUpdator.cxx:426
ProjectionMatricesSet.h
Trk::KalmanUpdator::calculateFilterStep
std::unique_ptr< TrackParameters > calculateFilterStep(const TrackParameters &, const Amg::Vector2D &, const Amg::MatrixX &, const int, FitQualityOnSurface *&, bool) const
Common maths calculation code for addToState and removeFromState - Amg::Vector2D interface.
Definition: KalmanUpdator.cxx:434
Trk::phi
@ phi
Definition: ParamDefs.h:75
Trk::KalmanUpdator::m_outputlevel
int m_outputlevel
MsgStream output level cached.
Definition: KalmanUpdator.h:254
AthAlgTool
Definition: AthAlgTool.h:26
Trk::KalmanUpdator::logInputCov
void logInputCov(const Amg::MatrixX &, const Amg::VectorX &, const Amg::MatrixX &) const
internal structuring: common logfile output of the inputs
Definition: KalmanUpdator.cxx:855
Amg::chi2
double chi2(const T &precision, const U &residual, const int sign=1)
Definition: EventPrimitivesCovarianceHelpers.h:221
Trk::KalmanUpdator::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)
Definition: KalmanUpdator.cxx:104
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37