CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/TrackingTools/KalmanUpdators/src/KFUpdator.cc

Go to the documentation of this file.
00001 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00002 #include "TrackingTools/PatternTools/interface/MeasurementExtractor.h"
00003 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00004 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
00005 #include "DataFormats/TrackingRecHit/interface/KfComponentsHolder.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "DataFormats/Math/interface/invertPosDefMatrix.h"
00008 #include "DataFormats/Math/interface/ProjectMatrix.h"
00009 
00010 TrajectoryStateOnSurface KFUpdator::update(const TrajectoryStateOnSurface& tsos,
00011                                            const TransientTrackingRecHit& aRecHit) const {
00012     switch (aRecHit.dimension()) {
00013         case 1: return update<1>(tsos,aRecHit);
00014         case 2: return update<2>(tsos,aRecHit);
00015         case 3: return update<3>(tsos,aRecHit);
00016         case 4: return update<4>(tsos,aRecHit);
00017         case 5: return update<5>(tsos,aRecHit);
00018     }
00019     throw cms::Exception("Rec hit of invalid dimension (not 1,2,3,4,5)") <<
00020          "The value was " << aRecHit.dimension() << 
00021         ", type is " << typeid(aRecHit).name() << 
00022         ", persistent rechit type " << (aRecHit.hit() ? typeid(*aRecHit.hit()).name() : "NULL") << "\n";
00023 }
00024 
00025 #define NEW
00026 #ifdef NEW
00027 template <unsigned int D>
00028 TrajectoryStateOnSurface KFUpdator::update(const TrajectoryStateOnSurface& tsos,
00029                                            const TransientTrackingRecHit& aRecHit) const {
00030 
00031   typedef typename AlgebraicROOTObject<D,5>::Matrix MatD5;
00032   typedef typename AlgebraicROOTObject<5,D>::Matrix Mat5D;
00033   typedef typename AlgebraicROOTObject<D,D>::SymMatrix SMatDD;
00034   typedef typename AlgebraicROOTObject<D>::Vector VecD;
00035   double pzSign = tsos.localParameters().pzSign();
00036 
00037   //MeasurementExtractor me(tsos);
00038 
00039   AlgebraicVector5 x = tsos.localParameters().vector();
00040   const AlgebraicSymMatrix55 &C = (tsos.localError().matrix());
00041   // Measurement matrix
00042   ProjectMatrix<double,5,D>  pf;
00043   MatD5 H; 
00044   VecD r, rMeas; 
00045   SMatDD V, VMeas;
00046 
00047   KfComponentsHolder holder; 
00048   holder.template setup<D>(&r, &V, &H, &pf, &rMeas, &VMeas, x, C);
00049   aRecHit.getKfComponents(holder);
00050   
00051   //MatD5 H = asSMatrix<D,5>(aRecHit.projectionMatrix());
00052 
00053   // Residuals of aPredictedState w.r.t. aRecHit, 
00054   //VecD r = asSVector<D>(aRecHit.parameters()) - me.measuredParameters<D>(aRecHit);
00055   //r -= me.measuredParameters<D>(aRecHit);
00056   r -= rMeas;
00057 
00058   // and covariance matrix of residuals
00059   //SMatDD V = asSMatrix<D>(aRecHit.parametersError());
00060   //SMatDD R = V + me.measuredError<D>(aRecHit);
00061   SMatDD R = V + VMeas;
00062   bool ok = invertPosDefMatrix(R);
00063   // error check moved at the end
00064   //R.Invert();
00065 
00066   // Compute Kalman gain matrix
00067   Mat5D K;
00068   AlgebraicMatrix55 M = AlgebraicMatrixID();
00069   if (holder.useProjFunc() ) {
00070     K = C*pf.project(R);
00071     pf.projectAndSubtractFrom(M,K);
00072   }
00073   else {
00074     K = (C * ROOT::Math::Transpose(H)) * R;
00075     M -=  K * H;
00076   }
00077 
00078   // Compute local filtered state vector
00079   AlgebraicVector5 fsv = x + K * r;
00080   // Compute covariance matrix of local filtered state vector
00081   AlgebraicSymMatrix55 fse = ROOT::Math::Similarity(M, C) + ROOT::Math::Similarity(K, V);
00082 
00083 
00084   /*
00085   // expanded similariy
00086   AlgebraicSymMatrix55 fse; 
00087   ROOT::Math::AssignSym::Evaluate(fse, (M* C) * ( ROOT::Math::Transpose(M)));
00088   AlgebraicSymMatrix55 tmp;
00089   ROOT::Math::AssignSym::Evaluate(tmp, (K*V) * (ROOT::Math::Transpose(K)));
00090   fse +=  tmp;
00091   */
00092 
00093   if (ok) {
00094     return TrajectoryStateOnSurface( LocalTrajectoryParameters(fsv, pzSign),
00095                                      LocalTrajectoryError(fse), tsos.surface(),&(tsos.globalParameters().magneticField()), tsos.surfaceSide() );
00096   }else {
00097     edm::LogError("KFUpdator")<<" could not invert martix:\n"<< (V+VMeas);
00098     return TrajectoryStateOnSurface();
00099   }
00100 
00101 }
00102 #endif
00103 
00104 #ifndef NEW
00105 template <unsigned int D>
00106 TrajectoryStateOnSurface KFUpdator::update(const TrajectoryStateOnSurface& tsos,
00107                                            const TransientTrackingRecHit& aRecHit) const {
00108 
00109   typedef typename AlgebraicROOTObject<D,5>::Matrix MatD5;
00110   typedef typename AlgebraicROOTObject<5,D>::Matrix Mat5D;
00111   typedef typename AlgebraicROOTObject<D,D>::SymMatrix SMatDD;
00112   typedef typename AlgebraicROOTObject<D>::Vector VecD;
00113   double pzSign = tsos.localParameters().pzSign();
00114 
00115   MeasurementExtractor me(tsos);
00116 
00117   AlgebraicVector5 x = tsos.localParameters().vector();
00118   const AlgebraicSymMatrix55 &C = (tsos.localError().matrix());
00119   // Measurement matrix
00120   MatD5 H = asSMatrix<D,5>(aRecHit.projectionMatrix());
00121 
00122   // Residuals of aPredictedState w.r.t. aRecHit, 
00123   VecD r = asSVector<D>(aRecHit.parameters()) - me.measuredParameters<D>(aRecHit);
00124 
00125   // and covariance matrix of residuals
00126   SMatDD V = asSMatrix<D>(aRecHit.parametersError());
00127   SMatDD R = V + me.measuredError<D>(aRecHit);
00128   int ierr = !  invertPosDefMatrix(R);; // if (ierr != 0) throw exception;
00129   //R.Invert();
00130 
00131   // Compute Kalman gain matrix
00132   Mat5D K = C * ROOT::Math::Transpose(H) * R;
00133 
00134   // Compute local filtered state vector
00135   AlgebraicVector5 fsv = x + K * r;
00136 
00137   // Compute covariance matrix of local filtered state vector
00138   AlgebraicMatrix55 I = AlgebraicMatrixID();
00139   AlgebraicMatrix55 M = I - K * H;
00140   AlgebraicSymMatrix55 fse = ROOT::Math::Similarity(M, C) + ROOT::Math::Similarity(K, V);
00141 
00142   return TrajectoryStateOnSurface( LocalTrajectoryParameters(fsv, pzSign),
00143                                    LocalTrajectoryError(fse), tsos.surface(),&(tsos.globalParameters().magneticField()), tsos.surfaceSide() );
00144 }
00145 #endif