CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KFUpdator.cc
Go to the documentation of this file.
9 
11  const TransientTrackingRecHit& aRecHit) const {
12  switch (aRecHit.dimension()) {
13  case 1: return update<1>(tsos,aRecHit);
14  case 2: return update<2>(tsos,aRecHit);
15  case 3: return update<3>(tsos,aRecHit);
16  case 4: return update<4>(tsos,aRecHit);
17  case 5: return update<5>(tsos,aRecHit);
18  }
19  throw cms::Exception("Rec hit of invalid dimension (not 1,2,3,4,5)") <<
20  "The value was " << aRecHit.dimension() <<
21  ", type is " << typeid(aRecHit).name() <<
22  ", persistent rechit type " << (aRecHit.hit() ? typeid(*aRecHit.hit()).name() : "NULL") << "\n";
23 }
24 
25 #define NEW
26 #ifdef NEW
27 template <unsigned int D>
29  const TransientTrackingRecHit& aRecHit) const {
30 
31  typedef typename AlgebraicROOTObject<D,5>::Matrix MatD5;
32  typedef typename AlgebraicROOTObject<5,D>::Matrix Mat5D;
33  typedef typename AlgebraicROOTObject<D,D>::SymMatrix SMatDD;
34  typedef typename AlgebraicROOTObject<D>::Vector VecD;
35  double pzSign = tsos.localParameters().pzSign();
36 
37  //MeasurementExtractor me(tsos);
38 
40  const AlgebraicSymMatrix55 &C = (tsos.localError().matrix());
41  // Measurement matrix
43  MatD5 H;
44  VecD r, rMeas;
45  SMatDD V, VMeas;
46 
47  KfComponentsHolder holder;
48  holder.template setup<D>(&r, &V, &H, &pf, &rMeas, &VMeas, x, C);
49  aRecHit.getKfComponents(holder);
50 
51  //MatD5 H = asSMatrix<D,5>(aRecHit.projectionMatrix());
52 
53  // Residuals of aPredictedState w.r.t. aRecHit,
54  //VecD r = asSVector<D>(aRecHit.parameters()) - me.measuredParameters<D>(aRecHit);
55  //r -= me.measuredParameters<D>(aRecHit);
56  r -= rMeas;
57 
58  // and covariance matrix of residuals
59  //SMatDD V = asSMatrix<D>(aRecHit.parametersError());
60  //SMatDD R = V + me.measuredError<D>(aRecHit);
61  SMatDD R = V + VMeas;
62  bool ok = invertPosDefMatrix(R);
63  // error check moved at the end
64  //R.Invert();
65 
66  // Compute Kalman gain matrix
67  Mat5D K;
69  if (holder.useProjFunc() ) {
70  K = C*pf.project(R);
71  pf.projectAndSubtractFrom(M,K);
72  }
73  else {
74  K = (C * ROOT::Math::Transpose(H)) * R;
75  M -= K * H;
76  }
77 
78  // Compute local filtered state vector
79  AlgebraicVector5 fsv = x + K * r;
80  // Compute covariance matrix of local filtered state vector
81  AlgebraicSymMatrix55 fse = ROOT::Math::Similarity(M, C) + ROOT::Math::Similarity(K, V);
82 
83 
84  /*
85  // expanded similariy
86  AlgebraicSymMatrix55 fse;
87  ROOT::Math::AssignSym::Evaluate(fse, (M* C) * ( ROOT::Math::Transpose(M)));
88  AlgebraicSymMatrix55 tmp;
89  ROOT::Math::AssignSym::Evaluate(tmp, (K*V) * (ROOT::Math::Transpose(K)));
90  fse += tmp;
91  */
92 
93  if (ok) {
95  LocalTrajectoryError(fse), tsos.surface(),&(tsos.globalParameters().magneticField()), tsos.surfaceSide() );
96  }else {
97  edm::LogError("KFUpdator")<<" could not invert martix:\n"<< (V+VMeas);
98  return TrajectoryStateOnSurface();
99  }
100 
101 }
102 #endif
103 
104 #ifndef NEW
105 template <unsigned int D>
107  const TransientTrackingRecHit& aRecHit) const {
108 
109  typedef typename AlgebraicROOTObject<D,5>::Matrix MatD5;
110  typedef typename AlgebraicROOTObject<5,D>::Matrix Mat5D;
111  typedef typename AlgebraicROOTObject<D,D>::SymMatrix SMatDD;
112  typedef typename AlgebraicROOTObject<D>::Vector VecD;
113  double pzSign = tsos.localParameters().pzSign();
114 
115  MeasurementExtractor me(tsos);
116 
118  const AlgebraicSymMatrix55 &C = (tsos.localError().matrix());
119  // Measurement matrix
120  MatD5 H = asSMatrix<D,5>(aRecHit.projectionMatrix());
121 
122  // Residuals of aPredictedState w.r.t. aRecHit,
123  VecD r = asSVector<D>(aRecHit.parameters()) - me.measuredParameters<D>(aRecHit);
124 
125  // and covariance matrix of residuals
126  SMatDD V = asSMatrix<D>(aRecHit.parametersError());
127  SMatDD R = V + me.measuredError<D>(aRecHit);
128  int ierr = ! invertPosDefMatrix(R);; // if (ierr != 0) throw exception;
129  //R.Invert();
130 
131  // Compute Kalman gain matrix
132  Mat5D K = C * ROOT::Math::Transpose(H) * R;
133 
134  // Compute local filtered state vector
135  AlgebraicVector5 fsv = x + K * r;
136 
137  // Compute covariance matrix of local filtered state vector
139  AlgebraicMatrix55 M = I - K * H;
140  AlgebraicSymMatrix55 fse = ROOT::Math::Similarity(M, C) + ROOT::Math::Similarity(K, V);
141 
143  LocalTrajectoryError(fse), tsos.surface(),&(tsos.globalParameters().magneticField()), tsos.surfaceSide() );
144 }
145 #endif
virtual int dimension() const =0
virtual void getKfComponents(KfComponentsHolder &holder) const
double pzSign() const
Sign of the z-component of the momentum in the local frame.
ROOT::Math::SMatrix< double, D1, D1, ROOT::Math::MatRepSym< double, D1 > > SymMatrix
bool useProjFunc() const
const LocalTrajectoryParameters & localParameters() const
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
virtual const TrackingRecHit * hit() const =0
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
virtual AlgebraicVector parameters() const =0
AlgebraicVector5 vector() const
ROOT::Math::SMatrix< double, D1, D2, ROOT::Math::MatRepStd< double, D1, D2 > > Matrix
const SurfaceType & surface() const
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
bool invertPosDefMatrix(ROOT::Math::SMatrix< T, N, N, ROOT::Math::MatRepSym< T, N > > &m)
const std::complex< double > I
Definition: I.h:8
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
ROOT::Math::SVector< double, D1 > Vector
ROOT::Math::SVector< double, 5 > AlgebraicVector5
const GlobalTrajectoryParameters & globalParameters() const
virtual AlgebraicMatrix projectionMatrix() const =0
#define update(a, b)
const MagneticField & magneticField() const
Definition: DDAxes.h:10
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:150
virtual AlgebraicSymMatrix parametersError() const =0