CMS 3D CMS Logo

CollinearFitAtTM.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFTracking/interface/CollinearFitAtTM.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 CollinearFitAtTM::CollinearFitAtTM ()
00005 {
00006   //
00007   // Jacobian
00008   //
00009   for ( int i=0; i<12; ++i ) {
00010     for ( int j=0; j<6; ++j )  jacobian_(i,j) = 0;
00011   }
00012   for ( int i=1; i<5; ++i ) {
00013     jacobian_(i,ParQpOut+i) = jacobian_(i+5,ParQpOut+i) = 1;
00014   }
00015   jacobian_(0,ParQpIn) = 1.;
00016   jacobian_(5,ParQpOut) = 1.;
00017 }
00018 
00019 bool
00020 CollinearFitAtTM::fit (const TrajectoryMeasurement& tm,
00021                        ResultVector& parameters,
00022                        ResultMatrix& covariance,
00023                        double& chi2)
00024 {
00025   //
00026   // check input
00027   //
00028   if ( !tm.forwardPredictedState().isValid() ||
00029        !tm.backwardPredictedState().isValid() ||
00030        !tm.updatedState().isValid() ) {
00031     edm::LogWarning("CollinearFitAtTM") << "Invalid state in TrajectoryMeasurement";
00032     return false;
00033   }
00034   //
00035   // prepare fit
00036   //
00037   AlgebraicVector5 fwdPar = tm.forwardPredictedState().localParameters().vector();
00038   AlgebraicSymMatrix55 fwdCov = tm.forwardPredictedState().localError().matrix();
00039   AlgebraicVector5 bwdPar = tm.backwardPredictedState().localParameters().vector();
00040   AlgebraicSymMatrix55 bwdCov = tm.backwardPredictedState().localError().matrix();
00041 
00042   LocalPoint hitPos(0.,0.,0.);
00043   LocalError hitErr(-1.,-1.,-1.);
00044   if ( tm.recHit()->isValid() ) {
00045     hitPos = tm.recHit()->localPosition();
00046     hitErr = tm.recHit()->localPositionError();
00047   }
00048 
00049   return fit(fwdPar,fwdCov,bwdPar,bwdCov,
00050              hitPos,hitErr,
00051              parameters,covariance,chi2);
00052 }
00053 
00054 bool
00055 CollinearFitAtTM::fit (const AlgebraicVector5& fwdParameters, 
00056                        const AlgebraicSymMatrix55& fwdCovariance,
00057                        const AlgebraicVector5& bwdParameters, 
00058                        const AlgebraicSymMatrix55& bwdCovariance,
00059                        const LocalPoint& hitPos, const LocalError& hitErr,
00060                        ResultVector& parameters, ResultMatrix& covariance,
00061                        double& chi2)
00062 {
00063 
00064   if ( hitErr.xx()>0 )
00065     jacobian_(10,ParX) = jacobian_(11,ParY) = 1;
00066   else
00067     jacobian_(10,ParX) = jacobian_(11,ParY) = 0;
00068 
00069   for ( int i=0; i<12; ++i ) {
00070     for ( int j=0; j<12; ++j )  weightMatrix_(i,j) = 0;
00071   }
00072 
00073   for ( int i=0; i<5; ++i )  measurements_(i) = fwdParameters(i);
00074   weightMatrix_.Place_at(fwdCovariance,0,0);
00075   for ( int i=0; i<5; ++i )  measurements_(i+5) = bwdParameters(i);
00076   weightMatrix_.Place_at(bwdCovariance,5,5);
00077   if ( hitErr.xx()>0 ) {
00078     measurements_(10) = hitPos.x();
00079     measurements_(11) = hitPos.y();
00080     weightMatrix_(10,10) = hitErr.xx();
00081     weightMatrix_(10,11) = weightMatrix_(11,10) = hitErr.xy();
00082     weightMatrix_(11,11) = hitErr.yy();
00083   }
00084   else {
00085     measurements_(10) = measurements_(11) = 0.;
00086     weightMatrix_(10,10) = weightMatrix_(11,11) = 1.;
00087     weightMatrix_(10,11) = weightMatrix_(11,10) = 0.;
00088   }
00089   //
00090   // invert covariance matrix
00091   //
00092   if ( !weightMatrix_.Invert() ) {
00093     edm::LogWarning("CollinearFitAtTM") << "Inversion of input covariance matrix failed";
00094     return false;
00095   }
00096 
00097   //
00098   projectedMeasurements_ = ROOT::Math::Transpose(jacobian_)*(weightMatrix_*measurements_);
00099   //
00100   // Fitted parameters and covariance matrix
00101   // 
00102   covariance = ROOT::Math::SimilarityT(jacobian_,weightMatrix_);
00103   if ( !covariance.Invert() ) {
00104     edm::LogWarning("CollinearFitAtTM") << "Inversion of resulting weight matrix failed";
00105     return false;
00106   }
00107 
00108   parameters = covariance*projectedMeasurements_;
00109 
00110   //
00111   // chi2
00112   //
00113   chi2 = ROOT::Math::Similarity(measurements_,weightMatrix_) -
00114     ROOT::Math::Similarity(projectedMeasurements_,covariance);
00115 
00116   return true;
00117 }

Generated on Tue Jun 9 17:44:49 2009 for CMSSW by  doxygen 1.5.4