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