CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TwoBodyDecayEstimator.cc
Go to the documentation of this file.
1 
4 
8 //#include "DataFormats/CLHEP/interface/Migration.h"
9 
11  :theNdf(0)
12 {
13  const edm::ParameterSet & estimatorConfig = config.getParameter< edm::ParameterSet >( "EstimatorParameters" );
14 
15  theRobustificationConstant = estimatorConfig.getUntrackedParameter< double >( "RobustificationConstant", 1.0 );
16  theMaxIterDiff = estimatorConfig.getUntrackedParameter< double >( "MaxIterationDifference", 1e-2 );
17  theMaxIterations = estimatorConfig.getUntrackedParameter< int >( "MaxIterations", 100 );
18  theUseInvariantMass = estimatorConfig.getUntrackedParameter< bool >( "UseInvariantMass", true );
19 }
20 
21 
22 TwoBodyDecay TwoBodyDecayEstimator::estimate( const std::vector< RefCountedLinearizedTrackState > & linTracks,
23  const TwoBodyDecayParameters & linearizationPoint,
24  const TwoBodyDecayVirtualMeasurement & vm ) const
25 {
26  if ( linTracks.size() != 2 )
27  {
28  edm::LogInfo( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::estimate"
29  << "Need 2 linearized tracks, got " << linTracks.size() << ".\n";
30  return TwoBodyDecay();
31  }
32 
33  AlgebraicVector vecM;
34  AlgebraicSymMatrix matG;
35  AlgebraicMatrix matA;
36 
37  bool check = constructMatrices( linTracks, linearizationPoint, vm, vecM, matG, matA );
38  if ( !check ) return TwoBodyDecay();
39 
40  AlgebraicSymMatrix matGPrime;
41  AlgebraicSymMatrix invAtGPrimeA;
42  AlgebraicVector vecEstimate;
43  AlgebraicVector res;
44 
45  int nIterations = 0;
46  bool stopIteration = false;
47 
48  // initialization - values are never used
49  int checkInversion = 0;
50  double chi2 = 0.;
51  double oldChi2 = 0.;
52  bool isValid = true;
53 
54  while( !stopIteration )
55  {
56  matGPrime = matG;
57 
58  // compute weights
59  if ( nIterations > 0 )
60  {
61  for ( int i = 0; i < 10; i++ )
62  {
63  double sigma = 1./sqrt( matG[i][i] );
64  double sigmaTimesR = sigma*theRobustificationConstant;
65  double absRes = fabs( res[i] );
66  if ( absRes > sigmaTimesR ) matGPrime[i][i] *= sigmaTimesR/absRes;
67  }
68  }
69 
70  // make LS-fit
71  invAtGPrimeA = ( matGPrime.similarityT(matA) ).inverse( checkInversion );
72  if ( checkInversion != 0 )
73  {
74  LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::estimate"
75  << "Matrix At*G'*A not invertible (in iteration " << nIterations
76  << ", ifail = " << checkInversion << ").\n";
77  isValid = false;
78  break;
79  }
80  vecEstimate = invAtGPrimeA*matA.T()*matGPrime*vecM;
81  res = matA*vecEstimate - vecM;
82  chi2 = dot( res, matGPrime*res );
83 
84  if ( ( nIterations > 0 ) && ( fabs( chi2 - oldChi2 ) < theMaxIterDiff ) ) stopIteration = true;
85  if ( nIterations == theMaxIterations ) stopIteration = true;
86 
87  oldChi2 = chi2;
88  nIterations++;
89  }
90 
91  if ( isValid )
92  {
93  AlgebraicSymMatrix pullsCov = matGPrime.inverse( checkInversion ) - invAtGPrimeA.similarity( matA );
94  thePulls = AlgebraicVector( matG.num_col(), 0 );
95  for ( int i = 0; i < pullsCov.num_col(); i++ ) thePulls[i] = res[i]/sqrt( pullsCov[i][i] );
96  }
97 
98  theNdf = matA.num_row() - matA.num_col();
99 
100  return TwoBodyDecay( TwoBodyDecayParameters( vecEstimate, invAtGPrimeA ), chi2, isValid, vm );
101 }
102 
103 
104 bool TwoBodyDecayEstimator::constructMatrices( const std::vector< RefCountedLinearizedTrackState > & linTracks,
105  const TwoBodyDecayParameters & linearizationPoint,
107  AlgebraicVector & vecM, AlgebraicSymMatrix & matG, AlgebraicMatrix & matA ) const
108 {
109 
110  PerigeeLinearizedTrackState* linTrack1 = dynamic_cast<PerigeeLinearizedTrackState*>( linTracks[0].get() );
111  PerigeeLinearizedTrackState* linTrack2 = dynamic_cast<PerigeeLinearizedTrackState*>( linTracks[1].get() );
112 
113  if (!linTrack1 || !linTrack2) return false;
114 
115  AlgebraicVector trackParam1 = asHepVector( linTrack1->predictedStateParameters() );
116  AlgebraicVector trackParam2 = asHepVector( linTrack2->predictedStateParameters() );
117 
118  if ( checkValues( trackParam1 ) || checkValues( trackParam2 ) || checkValues( linearizationPoint.parameters() ) ) return false;
119 
120  AlgebraicVector vecLinParam = linearizationPoint.sub( TwoBodyDecayParameters::px,
122 
123  double zMagField = linTrack1->track().field()->inInverseGeV( linTrack1->linearizationPoint() ).z();
124 
125  int checkInversion = 0;
126 
127  TwoBodyDecayDerivatives tpeDerivatives( linearizationPoint[TwoBodyDecayParameters::mass], vm.secondaryMass() );
128  std::pair< AlgebraicMatrix, AlgebraicMatrix > derivatives = tpeDerivatives.derivatives( linearizationPoint );
129 
130  TwoBodyDecayModel decayModel( linearizationPoint[TwoBodyDecayParameters::mass], vm.secondaryMass() );
131  std::pair< AlgebraicVector, AlgebraicVector > linCartMomenta = decayModel.cartesianSecondaryMomenta( linearizationPoint );
132 
133  // first track
134  AlgebraicMatrix matA1 = asHepMatrix( linTrack1->positionJacobian() );
135  AlgebraicMatrix matB1 = asHepMatrix( linTrack1->momentumJacobian() );
136  AlgebraicVector vecC1 = asHepVector( linTrack1->constantTerm() );
137 
138  AlgebraicVector curvMomentum1 = asHepVector( linTrack1->predictedStateMomentumParameters() );
139  AlgebraicMatrix curv2cart1 = decayModel.curvilinearToCartesianJacobian( curvMomentum1, zMagField );
140 
141  AlgebraicVector cartMomentum1 = decayModel.convertCurvilinearToCartesian( curvMomentum1, zMagField );
142  vecC1 += matB1*( curvMomentum1 - curv2cart1*cartMomentum1 );
143  matB1 = matB1*curv2cart1;
144 
145  AlgebraicMatrix matF1 = derivatives.first;
146  AlgebraicVector vecD1 = linCartMomenta.first - matF1*vecLinParam;
147  AlgebraicVector vecM1 = trackParam1 - vecC1 - matB1*vecD1;
148  AlgebraicSymMatrix covM1 = asHepMatrix( linTrack1->predictedStateError() );
149 
150 
151  AlgebraicSymMatrix matG1 = covM1.inverse( checkInversion );
152  if ( checkInversion != 0 )
153  {
154  LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::constructMatrices"
155  << "Matrix covM1 not invertible.";
156  return false;
157  }
158 
159  // diagonalize for robustification
160  AlgebraicMatrix matU1 = diagonalize( &matG1 ).T();
161 
162  // second track
163  AlgebraicMatrix matA2 = asHepMatrix( linTrack2->positionJacobian() );
164  AlgebraicMatrix matB2 = asHepMatrix( linTrack2->momentumJacobian() );
165  AlgebraicVector vecC2 = asHepVector( linTrack2->constantTerm() );
166 
167  AlgebraicVector curvMomentum2 = asHepVector( linTrack2->predictedStateMomentumParameters() );
168  AlgebraicMatrix curv2cart2 = decayModel.curvilinearToCartesianJacobian( curvMomentum2, zMagField );
169 
170  AlgebraicVector cartMomentum2 = decayModel.convertCurvilinearToCartesian( curvMomentum2, zMagField );
171  vecC2 += matB2*( curvMomentum2 - curv2cart2*cartMomentum2 );
172  matB2 = matB2*curv2cart2;
173 
174  AlgebraicMatrix matF2 = derivatives.second;
175  AlgebraicVector vecD2 = linCartMomenta.second - matF2*vecLinParam;
176  AlgebraicVector vecM2 = trackParam2 - vecC2 - matB2*vecD2;
177  AlgebraicSymMatrix covM2 = asHepMatrix( linTrack2->predictedStateError() );
178 
179  AlgebraicSymMatrix matG2 = covM2.inverse( checkInversion );
180  if ( checkInversion != 0 )
181  {
182  LogDebug( "Alignment" ) << "@SUB=TwoBodyDecayEstimator::constructMatrices"
183  << "Matrix covM2 not invertible.";
184  return false;
185  }
186 
187  // diagonalize for robustification
188  AlgebraicMatrix matU2 = diagonalize( &matG2 ).T();
189 
190  // combine first and second track
191  vecM = AlgebraicVector( 14, 0 );
192  vecM.sub( 1, matU1*vecM1 );
193  vecM.sub( 6, matU2*vecM2 );
194  // virtual measurement of the primary mass
195  vecM( 11 ) = vm.primaryMass();
196  // virtual measurement of the beam spot
197  vecM.sub( 12, vm.beamSpotPosition() );
198 
199  // full weight matrix
200  matG = AlgebraicSymMatrix( 14, 0 );
201  matG.sub( 1, matG1 );
202  matG.sub( 6, matG2 );
203  // virtual measurement error of the primary mass
204  matG[10][10] = 1./( vm.primaryWidth()*vm.primaryWidth() );
205  // virtual measurement error of the beam spot
206  matG.sub( 12, vm.beamSpotError().inverse( checkInversion ) );
207 
208  // full derivative matrix
209  matA = AlgebraicMatrix( 14, 9, 0 );
210  matA.sub( 1, 1, matU1*matA1 );
211  matA.sub( 6, 1, matU2*matA2 );
212  matA.sub( 1, 4, matU1*matB1*matF1 );
213  matA.sub( 6, 4, matU2*matB2*matF2 );
214  matA( 11, 9 ) = 1.;//mass
215  matA( 12, 1 ) = 1.;//vx
216  matA( 13, 2 ) = 1.;//vy
217  matA( 14, 3 ) = 1.;//vz
218 
219  return true;
220 }
221 
222 
224 {
225  bool isNan = false;
226  bool isInf = false;
227 
228  for ( int i = 0; i < vec.num_col(); ++i )
229  {
230  isNan = isNan || std::isnan( vec[i] );
231  isInf = isInf || std::isinf( vec[i] );
232  }
233 
234  return ( isNan || isInf );
235 }
#define LogDebug(id)
virtual AlgebraicVector3 predictedStateMomentumParameters() const
bool checkValues(const AlgebraicVector &vec) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
CLHEP::HepMatrix asHepMatrix(const ROOT::Math::SMatrix< double, N1, N2, typename ROOT::Math::MatRepStd< double, N1, N2 > > &rm)
Definition: Migration.h:49
const AlgebraicMatrix53 & momentumJacobian() const
const AlgebraicVector & parameters(void) const
Get decay parameters.
virtual reco::TransientTrack track() const
virtual TwoBodyDecay estimate(const std::vector< RefCountedLinearizedTrackState > &linTracks, const TwoBodyDecayParameters &linearizationPoint, const TwoBodyDecayVirtualMeasurement &vm) const
const MagneticField * field() const
double double double z
const std::pair< AlgebraicMatrix, AlgebraicMatrix > derivatives(const TwoBodyDecay &tbd) const
TwoBodyDecayEstimator(const edm::ParameterSet &config)
const AlgebraicVector beamSpotPosition(void) const
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:40
CLHEP::HepMatrix AlgebraicMatrix
bool isnan(float x)
Definition: math.h:13
T sqrt(T t)
Definition: SSEVec.h:46
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
const AlgebraicMatrix53 & positionJacobian() const
AlgebraicSymMatrix55 predictedStateError() const
const GlobalPoint & linearizationPoint() const
const AlgebraicSymMatrix beamSpotError(void) const
CLHEP::HepVector AlgebraicVector
const AlgebraicVector sub(ParameterName first, ParameterName last) const
Get specified range of decay parameters.
virtual bool constructMatrices(const std::vector< RefCountedLinearizedTrackState > &linTracks, const TwoBodyDecayParameters &linearizationPoint, const TwoBodyDecayVirtualMeasurement &vm, AlgebraicVector &vecM, AlgebraicSymMatrix &matG, AlgebraicMatrix &matA) const
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
CLHEP::HepSymMatrix AlgebraicSymMatrix
AlgebraicVector5 predictedStateParameters() const
CLHEP::HepVector asHepVector(const ROOT::Math::SVector< double, N > &v)
Definition: Migration.h:43
const std::pair< AlgebraicVector, AlgebraicVector > cartesianSecondaryMomenta(const AlgebraicVector &param)
const AlgebraicVector5 & constantTerm() const