CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/TwoBodyDecay/src/TwoBodyDecayEstimator.cc

Go to the documentation of this file.
00001 
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "MagneticField/Engine/interface/MagneticField.h"
00004 
00005 #include "Alignment/TwoBodyDecay/interface/TwoBodyDecayEstimator.h"
00006 #include "Alignment/TwoBodyDecay/interface/TwoBodyDecayModel.h"
00007 #include "Alignment/TwoBodyDecay/interface/TwoBodyDecayDerivatives.h"
00008 //#include "DataFormats/CLHEP/interface/Migration.h"
00009 
00010 TwoBodyDecayEstimator::TwoBodyDecayEstimator( const edm::ParameterSet & config )
00011   :theNdf(0)
00012 {
00013   const edm::ParameterSet & estimatorConfig = config.getParameter< edm::ParameterSet >( "EstimatorParameters" );
00014 
00015   theRobustificationConstant = estimatorConfig.getUntrackedParameter< double >( "RobustificationConstant", 1.0 );
00016   theMaxIterDiff = estimatorConfig.getUntrackedParameter< double >( "MaxIterationDifference", 1e-2 );
00017   theMaxIterations = estimatorConfig.getUntrackedParameter< int >( "MaxIterations", 100 );
00018   theUseInvariantMass = estimatorConfig.getUntrackedParameter< bool >( "UseInvariantMass", true );
00019 }
00020 
00021 
00022 TwoBodyDecay TwoBodyDecayEstimator::estimate( const std::vector< RefCountedLinearizedTrackState > & linTracks,
00023                                               const TwoBodyDecayParameters & linearizationPoint,
00024                                               const TwoBodyDecayVirtualMeasurement & vm ) const
00025 {
00026   if ( linTracks.size() != 2 )
00027   {
00028     edm::LogInfo( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::estimate"
00029                                 << "Need 2 linearized tracks, got " << linTracks.size() << ".\n";
00030     return TwoBodyDecay();
00031   }
00032 
00033   AlgebraicVector vecM;
00034   AlgebraicSymMatrix matG;
00035   AlgebraicMatrix matA;
00036 
00037   bool check = constructMatrices( linTracks, linearizationPoint, vm, vecM, matG, matA );
00038   if ( !check ) return TwoBodyDecay();
00039 
00040   AlgebraicSymMatrix matGPrime;
00041   AlgebraicSymMatrix invAtGPrimeA;
00042   AlgebraicVector vecEstimate;
00043   AlgebraicVector res;
00044 
00045   int nIterations = 0;
00046   bool stopIteration = false;
00047 
00048   // initialization - values are never used
00049   int checkInversion = 0;
00050   double chi2 = 0.;
00051   double oldChi2 = 0.;
00052   bool isValid = true;
00053 
00054   while( !stopIteration )
00055   {
00056     matGPrime = matG;
00057 
00058     // compute weights
00059     if ( nIterations > 0 )
00060     {
00061       for ( int i = 0; i < 10; i++ )
00062       {
00063         double sigma = 1./sqrt( matG[i][i] );
00064         double sigmaTimesR = sigma*theRobustificationConstant;
00065         double absRes = fabs( res[i] ); 
00066         if (  absRes > sigmaTimesR ) matGPrime[i][i] *= sigmaTimesR/absRes;
00067       }
00068     }
00069 
00070     // make LS-fit
00071     invAtGPrimeA = ( matGPrime.similarityT(matA) ).inverse( checkInversion );
00072     if ( checkInversion != 0 )
00073     {
00074       LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::estimate"
00075                               << "Matrix At*G'*A not invertible (in iteration " << nIterations
00076                               << ", ifail = " << checkInversion << ").\n";
00077       isValid = false;
00078       break;
00079     }
00080     vecEstimate = invAtGPrimeA*matA.T()*matGPrime*vecM;
00081     res = matA*vecEstimate - vecM;
00082     chi2 = dot( res, matGPrime*res );
00083 
00084     if ( ( nIterations > 0 ) && ( fabs( chi2 - oldChi2 ) < theMaxIterDiff ) ) stopIteration = true;
00085     if ( nIterations == theMaxIterations ) stopIteration = true;
00086 
00087     oldChi2 = chi2;
00088     nIterations++;
00089   }
00090 
00091   if ( isValid )
00092   {
00093     AlgebraicSymMatrix pullsCov = matGPrime.inverse( checkInversion ) - invAtGPrimeA.similarity( matA );
00094     thePulls = AlgebraicVector( matG.num_col(), 0 );
00095     for ( int i = 0; i < pullsCov.num_col(); i++ ) thePulls[i] = res[i]/sqrt( pullsCov[i][i] );
00096   }
00097 
00098   theNdf = matA.num_row() - matA.num_col();
00099 
00100   return TwoBodyDecay( TwoBodyDecayParameters( vecEstimate, invAtGPrimeA ), chi2, isValid, vm );
00101 }
00102 
00103 
00104 bool TwoBodyDecayEstimator::constructMatrices( const std::vector< RefCountedLinearizedTrackState > & linTracks,
00105                                                const TwoBodyDecayParameters & linearizationPoint,
00106                                                const TwoBodyDecayVirtualMeasurement & vm,
00107                                                AlgebraicVector & vecM, AlgebraicSymMatrix & matG, AlgebraicMatrix & matA ) const
00108 {
00109 
00110   PerigeeLinearizedTrackState* linTrack1 = dynamic_cast<PerigeeLinearizedTrackState*>( linTracks[0].get() );
00111   PerigeeLinearizedTrackState* linTrack2 = dynamic_cast<PerigeeLinearizedTrackState*>( linTracks[1].get() );
00112 
00113   if (!linTrack1 || !linTrack2) return false;
00114 
00115   AlgebraicVector trackParam1 = asHepVector( linTrack1->predictedStateParameters() );
00116   AlgebraicVector trackParam2 = asHepVector( linTrack2->predictedStateParameters() );
00117 
00118   if ( checkValues( trackParam1 ) || checkValues( trackParam2 ) || checkValues( linearizationPoint.parameters() ) ) return false;
00119 
00120   AlgebraicVector vecLinParam = linearizationPoint.sub( TwoBodyDecayParameters::px,
00121                                                         TwoBodyDecayParameters::mass );
00122 
00123   double zMagField = linTrack1->track().field()->inInverseGeV( linTrack1->linearizationPoint() ).z();
00124 
00125   int checkInversion = 0;
00126 
00127   TwoBodyDecayDerivatives tpeDerivatives( linearizationPoint[TwoBodyDecayParameters::mass], vm.secondaryMass() );
00128   std::pair< AlgebraicMatrix, AlgebraicMatrix > derivatives = tpeDerivatives.derivatives( linearizationPoint );
00129 
00130   TwoBodyDecayModel decayModel( linearizationPoint[TwoBodyDecayParameters::mass], vm.secondaryMass() );
00131   std::pair< AlgebraicVector, AlgebraicVector > linCartMomenta = decayModel.cartesianSecondaryMomenta( linearizationPoint );
00132 
00133   // first track
00134   AlgebraicMatrix matA1 = asHepMatrix( linTrack1->positionJacobian() );
00135   AlgebraicMatrix matB1 = asHepMatrix( linTrack1->momentumJacobian() );
00136   AlgebraicVector vecC1 = asHepVector( linTrack1->constantTerm() );
00137 
00138   AlgebraicVector curvMomentum1 = asHepVector( linTrack1->predictedStateMomentumParameters() );
00139   AlgebraicMatrix curv2cart1 = decayModel.curvilinearToCartesianJacobian( curvMomentum1, zMagField );
00140 
00141   AlgebraicVector cartMomentum1 = decayModel.convertCurvilinearToCartesian( curvMomentum1, zMagField );
00142   vecC1 += matB1*( curvMomentum1 - curv2cart1*cartMomentum1 );
00143   matB1 = matB1*curv2cart1;
00144 
00145   AlgebraicMatrix matF1 = derivatives.first;
00146   AlgebraicVector vecD1 = linCartMomenta.first - matF1*vecLinParam;
00147   AlgebraicVector vecM1 = trackParam1 - vecC1 - matB1*vecD1;
00148   AlgebraicSymMatrix covM1 = asHepMatrix( linTrack1->predictedStateError() );
00149 
00150 
00151   AlgebraicSymMatrix matG1 = covM1.inverse( checkInversion );
00152   if ( checkInversion != 0 )
00153   {
00154     LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::constructMatrices"
00155                             << "Matrix covM1 not invertible.";
00156     return false;
00157   }
00158 
00159   // diagonalize for robustification
00160    AlgebraicMatrix matU1 = diagonalize( &matG1 ).T();
00161 
00162   // second track
00163   AlgebraicMatrix matA2 = asHepMatrix( linTrack2->positionJacobian() );
00164   AlgebraicMatrix matB2 = asHepMatrix( linTrack2->momentumJacobian() );
00165   AlgebraicVector vecC2 = asHepVector( linTrack2->constantTerm() );
00166 
00167   AlgebraicVector curvMomentum2 = asHepVector( linTrack2->predictedStateMomentumParameters() );
00168   AlgebraicMatrix curv2cart2 = decayModel.curvilinearToCartesianJacobian( curvMomentum2, zMagField );
00169 
00170   AlgebraicVector cartMomentum2 = decayModel.convertCurvilinearToCartesian( curvMomentum2, zMagField );
00171   vecC2 += matB2*( curvMomentum2 - curv2cart2*cartMomentum2 );
00172   matB2 = matB2*curv2cart2;
00173 
00174   AlgebraicMatrix matF2 = derivatives.second;
00175   AlgebraicVector vecD2 = linCartMomenta.second - matF2*vecLinParam;
00176   AlgebraicVector vecM2 = trackParam2 - vecC2 - matB2*vecD2;
00177   AlgebraicSymMatrix covM2 = asHepMatrix( linTrack2->predictedStateError() );
00178 
00179   AlgebraicSymMatrix matG2 = covM2.inverse( checkInversion );
00180   if ( checkInversion != 0 )
00181   {
00182     LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::constructMatrices"
00183                             << "Matrix covM2 not invertible.";
00184     return false;
00185   }
00186 
00187   // diagonalize for robustification
00188   AlgebraicMatrix matU2 = diagonalize( &matG2 ).T();
00189 
00190   // combine first and second track
00191   vecM = AlgebraicVector( 14, 0 );
00192   vecM.sub( 1, matU1*vecM1 );
00193   vecM.sub( 6, matU2*vecM2 );
00194   // virtual measurement of the primary mass
00195   vecM( 11 ) = vm.primaryMass();
00196   // virtual measurement of the beam spot
00197   vecM.sub( 12, vm.beamSpotPosition() );
00198 
00199   // full weight matrix
00200   matG = AlgebraicSymMatrix( 14, 0 );
00201   matG.sub( 1, matG1 );
00202   matG.sub( 6, matG2 );
00203   // virtual measurement error of the primary mass
00204   matG[10][10] = 1./( vm.primaryWidth()*vm.primaryWidth() );
00205   // virtual measurement error of the beam spot
00206   matG.sub( 12, vm.beamSpotError().inverse( checkInversion ) );
00207 
00208   // full derivative matrix
00209   matA = AlgebraicMatrix( 14, 9, 0 );
00210   matA.sub( 1, 1, matU1*matA1 );
00211   matA.sub( 6, 1, matU2*matA2 );
00212   matA.sub( 1, 4, matU1*matB1*matF1 );
00213   matA.sub( 6, 4, matU2*matB2*matF2 );
00214   matA( 11, 9 ) = 1.;//mass
00215   matA( 12, 1 ) = 1.;//vx
00216   matA( 13, 2 ) = 1.;//vy
00217   matA( 14, 3 ) = 1.;//vz
00218 
00219   return true;
00220 }
00221 
00222 
00223 bool TwoBodyDecayEstimator::checkValues( const AlgebraicVector & vec ) const
00224 {
00225   bool isNan = false;
00226   bool isInf = false;
00227 
00228   for ( int i = 0; i < vec.num_col(); ++i )
00229   {
00230     isNan = isNan || std::isnan( vec[i] );
00231     isInf = isInf || std::isinf( vec[i] );
00232   }
00233 
00234   return ( isNan || isInf );
00235 }