CMS 3D CMS Logo

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