41 theLabel(0), theOffset(0), measDim(0), transFlag(false), measTransformation(), scatFlag(
42 false), localDerivatives(), globalLabels(), globalDerivatives() {
44 for (
unsigned int i = 0; i < 5; ++i) {
45 for (
unsigned int j = 0; j < 5; ++j) {
52 theLabel(0), theOffset(0), p2pJacobian(aJacobian), measDim(0), transFlag(
53 false), measTransformation(), scatFlag(false), localDerivatives(), globalLabels(), globalDerivatives() {
57 GblPoint::~GblPoint() {
70 const TVectorD &aResiduals,
const TVectorD &aPrecision,
71 double minPrecision) {
72 measDim = aResiduals.GetNrows();
73 unsigned int iOff = 5 -
measDim;
74 for (
unsigned int i = 0; i <
measDim; ++i) {
77 aPrecision[i] >= minPrecision ? aPrecision[i] : 0.);
78 for (
unsigned int j = 0; j <
measDim; ++j) {
95 const TVectorD &aResiduals,
const TMatrixDSym &aPrecision,
96 double minPrecision) {
97 measDim = aResiduals.GetNrows();
98 TMatrixDSymEigen measEigen(aPrecision);
104 TVectorD transPrecision = measEigen.GetEigenValues();
106 unsigned int iOff = 5 -
measDim;
107 for (
unsigned int i = 0; i <
measDim; ++i) {
110 transPrecision[i] >= minPrecision ? transPrecision[i] : 0.);
111 for (
unsigned int j = 0; j <
measDim; ++j) {
126 const TVectorD &aPrecision,
double minPrecision) {
127 measDim = aResiduals.GetNrows();
128 unsigned int iOff = 5 -
measDim;
129 for (
unsigned int i = 0; i <
measDim; ++i) {
132 aPrecision[i] >= minPrecision ? aPrecision[i] : 0.);
147 const TMatrixDSym &aPrecision,
double minPrecision) {
148 measDim = aResiduals.GetNrows();
149 TMatrixDSymEigen measEigen(aPrecision);
155 TVectorD transPrecision = measEigen.GetEigenValues();
156 unsigned int iOff = 5 -
measDim;
157 for (
unsigned int i = 0; i <
measDim; ++i) {
160 transPrecision[i] >= minPrecision ? transPrecision[i] : 0.);
161 for (
unsigned int j = 0; j <
measDim; ++j) {
183 SVector5 &aPrecision)
const {
198 aTransformation.UnitMatrix();
211 const TVectorD &aPrecision) {
237 const TMatrixDSym &aPrecision) {
239 TMatrixDSymEigen scatEigen(aPrecision);
240 TMatrixD aTransformation = scatEigen.GetEigenVectors();
242 TVectorD transResiduals = aTransformation * aResiduals;
243 TVectorD transPrecision = scatEigen.GetEigenValues();
244 for (
unsigned int i = 0; i < 2; ++i) {
247 for (
unsigned int j = 0; j < 2; ++j) {
265 SVector2 &aPrecision)
const {
276 aTransformation.ResizeTo(2, 2);
278 for (
unsigned int i = 0; i < 2; ++i) {
279 for (
unsigned int j = 0; j < 2; ++j) {
284 aTransformation.UnitMatrix();
321 const TMatrixD &aDerivatives) {
389 SMatrix23 CA = aJac.Sub<SMatrix23>(3, 0)
390 * aJac.Sub<SMatrix33>(0, 0).InverseFast(ifail);
391 SMatrix22 DCAB = aJac.Sub<SMatrix22>(3, 3) - CA * aJac.Sub<SMatrix32>(0, 3);
416 SVector2 &vecWd)
const {
420 if (aDirection < 1) {
430 if (!matW.InvertFast()) {
431 std::cout <<
" GblPoint::getDerivatives failed to invert matrix: "
432 << matW << std::endl;
434 <<
" Possible reason for singular matrix: multiple GblPoints at same arc-length"
436 throw std::overflow_error(
"Singular matrix inversion exception");
448 std::cout <<
" GblPoint";
450 std::cout <<
", label " <<
theLabel;
456 std::cout <<
", " <<
measDim <<
" measurements";
459 std::cout <<
", scatterer";
462 std::cout <<
", diagonalized";
466 <<
" local derivatives";
470 <<
" global derivatives";
472 std::cout << std::endl;
475 std::cout <<
" Measurement" << std::endl;
482 std::cout <<
" Scatterer" << std::endl;
487 std::cout <<
" Local Derivatives:" << std::endl;
491 std::cout <<
" Global Labels:";
492 for (
unsigned int i = 0; i <
globalLabels.size(); ++i) {
495 std::cout << std::endl;
496 std::cout <<
" Global Derivatives:" << std::endl;
499 std::cout <<
" Jacobian " << std::endl;
500 std::cout <<
" Point-to-point " << std::endl <<
p2pJacobian
503 std::cout <<
" To previous offset " << std::endl <<
prevJacobian
505 std::cout <<
" To next offset " << std::endl <<
nextJacobian
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.