CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
TwoBodyDecayEstimator Class Reference

#include <TwoBodyDecayEstimator.h>

Public Types

typedef
PerigeeLinearizedTrackState::RefCountedLinearizedTrackState 
RefCountedLinearizedTrackState
 

Public Member Functions

virtual TwoBodyDecayEstimatorclone (void) const
 
virtual TwoBodyDecay estimate (const std::vector< RefCountedLinearizedTrackState > &linTracks, const TwoBodyDecayParameters &linearizationPoint, const TwoBodyDecayVirtualMeasurement &vm) const
 
int ndf (void) const
 
const AlgebraicVectorpulls (void) const
 
 TwoBodyDecayEstimator (const edm::ParameterSet &config)
 
virtual ~TwoBodyDecayEstimator (void)
 

Protected Member Functions

virtual bool constructMatrices (const std::vector< RefCountedLinearizedTrackState > &linTracks, const TwoBodyDecayParameters &linearizationPoint, const TwoBodyDecayVirtualMeasurement &vm, AlgebraicVector &vecM, AlgebraicSymMatrix &matG, AlgebraicMatrix &matA) const
 

Private Member Functions

bool checkValues (const AlgebraicVector &vec) const
 

Private Attributes

int theMaxIterations
 
double theMaxIterDiff
 
int theNdf
 
AlgebraicVector thePulls
 
double theRobustificationConstant
 
bool theUseInvariantMass
 

Detailed Description

Class TwoBodyDecayEstimator estimates the decay parameters and the corresponding error matrix (see TwoBodyDecayParameter) from two linearized tracks, and an hypothesis for the primary particle's mass and width as well as the mass of the secondary particles. It utilizes a robust M-estimator.

/author Edmund Widl

Definition at line 20 of file TwoBodyDecayEstimator.h.

Member Typedef Documentation

Definition at line 25 of file TwoBodyDecayEstimator.h.

Constructor & Destructor Documentation

TwoBodyDecayEstimator::TwoBodyDecayEstimator ( const edm::ParameterSet config)

Definition at line 11 of file TwoBodyDecayEstimator.cc.

References alignCSCRings::e, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), theMaxIterations, theMaxIterDiff, theRobustificationConstant, theUseInvariantMass, and funct::true.

Referenced by clone().

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 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual TwoBodyDecayEstimator::~TwoBodyDecayEstimator ( void  )
inlinevirtual

Definition at line 28 of file TwoBodyDecayEstimator.h.

28 {}

Member Function Documentation

bool TwoBodyDecayEstimator::checkValues ( const AlgebraicVector vec) const
private

Definition at line 224 of file TwoBodyDecayEstimator.cc.

References i, and edm::isNotFinite().

Referenced by constructMatrices().

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 }
int i
Definition: DBlmapReader.cc:9
bool isNotFinite(T x)
Definition: isFinite.h:10
virtual TwoBodyDecayEstimator* TwoBodyDecayEstimator::clone ( void  ) const
inlinevirtual

Definition at line 37 of file TwoBodyDecayEstimator.h.

References TwoBodyDecayEstimator().

37 { return new TwoBodyDecayEstimator( *this ); }
TwoBodyDecayEstimator(const edm::ParameterSet &config)
bool TwoBodyDecayEstimator::constructMatrices ( const std::vector< RefCountedLinearizedTrackState > &  linTracks,
const TwoBodyDecayParameters linearizationPoint,
const TwoBodyDecayVirtualMeasurement vm,
AlgebraicVector vecM,
AlgebraicSymMatrix matG,
AlgebraicMatrix matA 
) const
protectedvirtual

Definition at line 105 of file TwoBodyDecayEstimator.cc.

References asHepMatrix(), asHepVector(), TwoBodyDecayVirtualMeasurement::beamSpotError(), TwoBodyDecayVirtualMeasurement::beamSpotPosition(), TwoBodyDecayModel::cartesianSecondaryMomenta(), checkValues(), PerigeeLinearizedTrackState::constantTerm(), TwoBodyDecayDerivatives::derivatives(), reco::TransientTrack::field(), MagneticField::inInverseGeV(), PerigeeLinearizedTrackState::linearizationPoint(), LogDebug, TwoBodyDecayParameters::mass, PerigeeLinearizedTrackState::momentumJacobian(), TwoBodyDecayParameters::parameters(), PerigeeLinearizedTrackState::positionJacobian(), PerigeeLinearizedTrackState::predictedStateError(), PerigeeLinearizedTrackState::predictedStateMomentumParameters(), PerigeeLinearizedTrackState::predictedStateParameters(), TwoBodyDecayVirtualMeasurement::primaryMass(), TwoBodyDecayVirtualMeasurement::primaryWidth(), TwoBodyDecayParameters::px, TwoBodyDecayVirtualMeasurement::secondaryMass(), TwoBodyDecayParameters::sub(), PerigeeLinearizedTrackState::track(), and z.

Referenced by estimate().

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 }
#define LogDebug(id)
virtual AlgebraicVector3 predictedStateMomentumParameters() const
bool checkValues(const AlgebraicVector &vec) const
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
const MagneticField * field() const
const std::pair< AlgebraicMatrix, AlgebraicMatrix > derivatives(const TwoBodyDecay &tbd) const
const AlgebraicVector beamSpotPosition(void) const
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:39
CLHEP::HepMatrix AlgebraicMatrix
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.
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
TwoBodyDecay TwoBodyDecayEstimator::estimate ( const std::vector< RefCountedLinearizedTrackState > &  linTracks,
const TwoBodyDecayParameters linearizationPoint,
const TwoBodyDecayVirtualMeasurement vm 
) const
virtual

Definition at line 23 of file TwoBodyDecayEstimator.cc.

References DDVectorGetter::check(), constructMatrices(), dot(), i, LogDebug, mathSSE::sqrt(), theMaxIterations, theMaxIterDiff, theNdf, thePulls, and theRobustificationConstant.

Referenced by TwoBodyDecayFitter::estimate().

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 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
CLHEP::HepMatrix AlgebraicMatrix
bool check(const std::string &)
T sqrt(T t)
Definition: SSEVec.h:48
CLHEP::HepVector AlgebraicVector
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
int TwoBodyDecayEstimator::ndf ( void  ) const
inline

Definition at line 34 of file TwoBodyDecayEstimator.h.

References theNdf.

34 { return theNdf; }
const AlgebraicVector& TwoBodyDecayEstimator::pulls ( void  ) const
inline

Definition at line 35 of file TwoBodyDecayEstimator.h.

References thePulls.

35 { return thePulls; }

Member Data Documentation

int TwoBodyDecayEstimator::theMaxIterations
private

Definition at line 52 of file TwoBodyDecayEstimator.h.

Referenced by estimate(), and TwoBodyDecayEstimator().

double TwoBodyDecayEstimator::theMaxIterDiff
private

Definition at line 51 of file TwoBodyDecayEstimator.h.

Referenced by estimate(), and TwoBodyDecayEstimator().

int TwoBodyDecayEstimator::theNdf
mutableprivate

Definition at line 55 of file TwoBodyDecayEstimator.h.

Referenced by estimate(), and ndf().

AlgebraicVector TwoBodyDecayEstimator::thePulls
mutableprivate

Definition at line 56 of file TwoBodyDecayEstimator.h.

Referenced by estimate(), and pulls().

double TwoBodyDecayEstimator::theRobustificationConstant
private

Definition at line 50 of file TwoBodyDecayEstimator.h.

Referenced by estimate(), and TwoBodyDecayEstimator().

bool TwoBodyDecayEstimator::theUseInvariantMass
private

Definition at line 53 of file TwoBodyDecayEstimator.h.

Referenced by TwoBodyDecayEstimator().