39 #include "TMatrixDSym.h"
40 #include "TMatrixDSymEigen.h"
42 #include "Math/SMatrix.h"
43 #include "Math/SVector.h"
44 typedef ROOT::Math::SMatrix<double, 2> SMatrix22;
45 typedef ROOT::Math::SMatrix<double, 2, 3> SMatrix23;
46 typedef ROOT::Math::SMatrix<double, 2, 5> SMatrix25;
47 typedef ROOT::Math::SMatrix<double, 2, 7> SMatrix27;
48 typedef ROOT::Math::SMatrix<double, 3, 2> SMatrix32;
49 typedef ROOT::Math::SMatrix<double, 3> SMatrix33;
50 typedef ROOT::Math::SMatrix<double, 5> SMatrix55;
51 typedef ROOT::Math::SVector<double, 2> SVector2;
52 typedef ROOT::Math::SVector<double, 5> SVector5;
70 explicit GblPoint(
const TMatrixD &aJacobian);
71 explicit GblPoint(
const SMatrix55 &aJacobian);
73 void addMeasurement(
const TMatrixD &aProjection,
const TVectorD &aResiduals,
74 const TVectorD &aPrecision,
double minPrecision = 0.);
75 void addMeasurement(
const TMatrixD &aProjection,
const TVectorD &aResiduals,
76 const TMatrixDSym &aPrecision,
double minPrecision = 0.);
77 void addMeasurement(
const TVectorD &aResiduals,
const TVectorD &aPrecision,
78 double minPrecision = 0.);
80 const TMatrixDSym &aPrecision,
double minPrecision = 0.);
83 SVector5 &aPrecision)
const;
85 void addScatterer(
const TVectorD &aResiduals,
const TVectorD &aPrecision);
87 const TMatrixDSym &aPrecision);
89 void getScatterer(SMatrix22 &aTransformation, SVector2 &aResiduals,
90 SVector2 &aPrecision)
const;
92 void addLocals(
const TMatrixD &aDerivatives);
95 void addGlobals(
const std::vector<int> &aLabels,
96 const TMatrixD &aDerivatives);
107 void getDerivatives(
int aDirection, SMatrix22 &matW, SMatrix22 &matWJ,
108 SVector2 &vecWd)
const;
109 void printPoint(
unsigned int level = 0)
const;
unsigned int getNumGlobals() const
Retrieve number of global derivatives from a point.
unsigned int hasMeasurement() const
Check for measurement at a point.
int getOffset() const
Retrieve offset for point.
void addMeasurement(const TMatrixD &aProjection, const TVectorD &aResiduals, const TVectorD &aPrecision, double minPrecision=0.)
Add a measurement to a point.
void addNextJacobian(const SMatrix55 &aJac)
Define jacobian to next scatterer (by GBLTrajectory constructor)
unsigned int measDim
Dimension of measurement (1-5), 0 indicates absence of measurement.
SMatrix55 measProjection
Projection from measurement to local system.
const TMatrixD & getLocalDerivatives() const
Retrieve local derivatives from a point.
SVector5 measResiduals
Measurement residuals.
SMatrix55 nextJacobian
Jacobian to next scatterer (or last measurement)
unsigned int getLabel() const
Retrieve label of point.
bool scatFlag
Scatterer present?
GblPoint(const TMatrixD &aJacobian)
Create a point.
int theOffset
Offset number at point if not negative (else interpolation needed)
SVector2 scatResiduals
Scattering residuals (initial kinks if iterating)
void getDerivatives(int aDirection, SMatrix22 &matW, SMatrix22 &matWJ, SVector2 &vecWd) const
Retrieve derivatives of local track model.
SVector5 measPrecision
Measurement precision (diagonal of inverse covariance matrix)
const TMatrixD & getGlobalDerivatives() const
Retrieve global derivatives from a point.
void getScatterer(SMatrix22 &aTransformation, SVector2 &aResiduals, SVector2 &aPrecision) const
Retrieve scatterer of a point.
SMatrix22 scatTransformation
Transformation of diagonalization (of scat. precision matrix)
unsigned int theLabel
Label identifying point.
SMatrix55 prevJacobian
Jacobian to previous scatterer (or first measurement)
void setOffset(int anOffset)
Define offset for point (by GBLTrajectory constructor)
SMatrix55 p2pJacobian
Point-to-point jacobian from previous point.
unsigned int getNumLocals() const
Retrieve number of local derivatives from a point.
std::vector< int > globalLabels
Labels of global (MP-II) derivatives.
SVector2 scatPrecision
Scattering precision (diagonal of inverse covariance matrix)
void getMeasTransformation(TMatrixD &aTransformation) const
Get measurement transformation (from diagonalization).
void addPrevJacobian(const SMatrix55 &aJac)
Define jacobian to previous scatterer (by GBLTrajectory constructor)
void printPoint(unsigned int level=0) const
Print GblPoint.
void getScatTransformation(TMatrixD &aTransformation) const
Get scatterer transformation (from diagonalization).
const SMatrix55 & getP2pJacobian() const
Retrieve point-to-(previous)point jacobian.
TMatrixD localDerivatives
Derivatives of measurement vs additional local (fit) parameters.
TMatrixD globalDerivatives
Derivatives of measurement vs additional global (MP-II) parameters.
void getMeasurement(SMatrix55 &aProjection, SVector5 &aResiduals, SVector5 &aPrecision) const
Retrieve measurement of a point.
bool transFlag
Transformation exists?
void addGlobals(const std::vector< int > &aLabels, const TMatrixD &aDerivatives)
Add global derivatives to a point.
void addScatterer(const TVectorD &aResiduals, const TVectorD &aPrecision)
Add a (thin) scatterer to a point.
bool hasScatterer() const
Check for scatterer at a point.
void addLocals(const TMatrixD &aDerivatives)
Add local derivatives to a point.
std::vector< int > getGlobalLabels() const
Retrieve global derivatives labels from a point.
void setLabel(unsigned int aLabel)
Define label of point (by GBLTrajectory constructor)
TMatrixD measTransformation
Transformation of diagonalization (of meas. precision matrix)
Namespace for the general broken lines package.