CMS 3D CMS Logo

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