GeneralBrokenLines  V02-02-01
using EIGEN
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
gbl::GblPoint Class Reference

Point on trajectory. More...

#include <GblPoint.h>

Public Member Functions

 GblPoint (const Matrix5d &aJacobian)
 Create a point. More...
 
 GblPoint (const GblPoint &)=default
 
GblPointoperator= (const GblPoint &)=default
 
 GblPoint (GblPoint &&)=default
 
GblPointoperator= (GblPoint &&)=default
 
virtual ~GblPoint ()
 
template<typename Projection , typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime !=1)>::type * = nullptr>
void addMeasurement (const Eigen::MatrixBase< Projection > &aProjection, const Eigen::MatrixBase< Residuals > &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision, double minPrecision=0.)
 Add a measurement to a point. More...
 
template<typename Projection , typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void addMeasurement (const Eigen::MatrixBase< Projection > &aProjection, const Eigen::MatrixBase< Residuals > &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision, double minPrecision=0.)
 Add a measurement to a point. More...
 
template<typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime !=1)>::type * = nullptr>
void addMeasurement (const Eigen::MatrixBase< Residuals > &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision, double minPrecision=0.)
 Add a measurement to a point. More...
 
template<typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void addMeasurement (const Eigen::MatrixBase< Residuals > &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision, double minPrecision=0.)
 Add a measurement to a point. More...
 
template<typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==2)>::type * = nullptr>
void addScatterer (const Eigen::Vector2d &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision)
 Add a (thin) scatterer to a point. More...
 
template<typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void addScatterer (const Eigen::Vector2d &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision)
 Add a (thin) scatterer to a point. More...
 
template<typename Derivative >
void addLocals (const Eigen::MatrixBase< Derivative > &aDerivatives)
 Add local derivatives to a point. More...
 
template<typename Derivative >
void addGlobals (const std::vector< int > &aLabels, const Eigen::MatrixBase< Derivative > &aDerivatives)
 Add global derivatives to a point. More...
 
unsigned int hasMeasurement () const
 Check for measurement at a point. More...
 
double getMeasPrecMin () const
 get precision cutoff. More...
 
void getMeasurement (Matrix5d &aProjection, Vector5d &aResiduals, Vector5d &aPrecision) const
 Retrieve measurement of a point. More...
 
void getMeasTransformation (Eigen::MatrixXd &aTransformation) const
 Get measurement transformation (from diagonalization). More...
 
bool hasScatterer () const
 Check for scatterer at a point. More...
 
void getScatterer (Eigen::Matrix2d &aTransformation, Eigen::Vector2d &aResiduals, Eigen::Vector2d &aPrecision) const
 Retrieve scatterer of a point. More...
 
void getScatTransformation (Eigen::Matrix2d &aTransformation) const
 Get scatterer transformation (from diagonalization). More...
 
unsigned int getNumLocals () const
 Retrieve number of local derivatives from a point. More...
 
const Eigen::MatrixXd & getLocalDerivatives () const
 Retrieve local derivatives from a point. More...
 
unsigned int getNumGlobals () const
 Retrieve number of global derivatives from a point. More...
 
void getGlobalLabels (std::vector< int > &aLabels) const
 Retrieve global derivatives labels from a point. More...
 
void getGlobalDerivatives (Eigen::MatrixXd &aDerivatives) const
 Retrieve global derivatives from a point. More...
 
void getGlobalLabelsAndDerivatives (unsigned int aRow, std::vector< int > &aLabels, std::vector< double > &aDerivatives) const
 Retrieve global derivatives from a point for a single row. More...
 
unsigned int getLabel () const
 Retrieve label of point. More...
 
int getOffset () const
 Retrieve offset for point. More...
 
const Matrix5dgetP2pJacobian () const
 Retrieve point-to-(previous)point jacobian. More...
 
void getDerivatives (int aDirection, Eigen::Matrix2d &matW, Eigen::Matrix2d &matWJ, Eigen::Vector2d &vecWd) const
 Retrieve derivatives of local track model. More...
 
void printPoint (unsigned int level=0) const
 Print GblPoint. More...
 

Private Member Functions

void setLabel (unsigned int aLabel)
 Define label of point (by GBLTrajectory constructor) More...
 
void setOffset (int anOffset)
 Define offset for point (by GBLTrajectory constructor) More...
 
void addPrevJacobian (const Matrix5d &aJac)
 Define jacobian to previous scatterer (by GBLTrajectory constructor) More...
 
void addNextJacobian (const Matrix5d &aJac)
 Define jacobian to next scatterer (by GBLTrajectory constructor) More...
 

Private Attributes

unsigned int theLabel
 Label identifying point. More...
 
int theOffset
 Offset number at point if not negative (else interpolation needed) More...
 
Matrix5d p2pJacobian
 Point-to-point jacobian from previous point. More...
 
Matrix5d prevJacobian
 Jacobian to previous scatterer (or first measurement) More...
 
Matrix5d nextJacobian
 Jacobian to next scatterer (or last measurement) More...
 
unsigned int measDim
 Dimension of measurement (1-5), 0 indicates absence of measurement. More...
 
double measPrecMin
 Minimal measurement precision (for usage) More...
 
Matrix5d measProjection
 Projection from measurement to local system. More...
 
Vector5d measResiduals
 Measurement residuals. More...
 
Vector5d measPrecision
 Measurement precision (diagonal of inverse covariance matrix) More...
 
bool transFlag
 Transformation exists? More...
 
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor, 5, 5 > measTransformation
 Transformation of diagonalization (of meas. precision matrix) More...
 
bool scatFlag
 Scatterer present? More...
 
Eigen::Matrix2d scatTransformation
 Transformation of diagonalization (of scat. precision matrix) More...
 
Eigen::Vector2d scatResiduals
 Scattering residuals (initial kinks if iterating) More...
 
Eigen::Vector2d scatPrecision
 Scattering precision (diagonal of inverse covariance matrix) More...
 
Eigen::MatrixXd localDerivatives
 Derivatives of measurement vs additional local (fit) parameters. More...
 
std::vector< int > globalLabels
 Labels of global (MP-II) derivatives. More...
 
Eigen::MatrixXd globalDerivatives
 Derivatives of measurement vs additional global (MP-II) parameters. More...
 

Friends

class GblTrajectory
 

Detailed Description

Point on trajectory.

User supplied point on (initial) trajectory.

Must have jacobian for propagation from previous point. May have:

  1. Measurement (1D - 5D)
  2. Scatterer (thin, 2D kinks)
  3. Additional local parameters (with derivatives). Fitted together with track parameters.
  4. Additional global parameters (with labels and derivatives). Not fitted, only passed on to (binary) file for fitting with Millepede-II.

Definition at line 69 of file GblPoint.h.

Constructor & Destructor Documentation

◆ GblPoint() [1/3]

gbl::GblPoint::GblPoint ( const Matrix5d aJacobian)

Create a point.

Create point on (initial) trajectory. Needs transformation jacobian from previous point.

Parameters
[in]aJacobianTransformation jacobian from previous point

Definition at line 41 of file GblPoint.cpp.

◆ GblPoint() [2/3]

gbl::GblPoint::GblPoint ( const GblPoint )
default

◆ GblPoint() [3/3]

gbl::GblPoint::GblPoint ( GblPoint &&  )
default

◆ ~GblPoint()

gbl::GblPoint::~GblPoint ( )
virtual

Definition at line 65 of file GblPoint.cpp.

Member Function Documentation

◆ addGlobals()

template<typename Derivative >
void gbl::GblPoint::addGlobals ( const std::vector< int > &  aLabels,
const Eigen::MatrixBase< Derivative > &  aDerivatives 
)

Add global derivatives to a point.

Point needs to have a measurement.

Template Parameters
DerivativeDerivatives matrix
Parameters
[in]aLabelsGlobal derivatives labels
[in]aDerivativesGlobal derivatives (matrix)

Definition at line 403 of file GblPoint.h.

References globalDerivatives, globalLabels, measDim, measTransformation, and transFlag.

Referenced by exampleDc(), and exampleSit().

◆ addLocals()

template<typename Derivative >
void gbl::GblPoint::addLocals ( const Eigen::MatrixBase< Derivative > &  aDerivatives)

Add local derivatives to a point.

Point needs to have a measurement.

Template Parameters
DerivativeDerivatives matrix
Parameters
[in]aDerivativesLocal derivatives (matrix)

Definition at line 391 of file GblPoint.h.

References localDerivatives, measDim, measTransformation, and transFlag.

Referenced by example3().

◆ addMeasurement() [1/4]

template<typename Projection , typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * >
void gbl::GblPoint::addMeasurement ( const Eigen::MatrixBase< Projection > &  aProjection,
const Eigen::MatrixBase< Residuals > &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision,
double  minPrecision = 0. 
)

Add a measurement to a point.

Add measurement (in measurement system) with arbitrary precision (inverse covariance) matrix. Will be diagonalized. ((up to) 2D: position, 4D: slope+position, 5D: curvature+slope+position)

Template Parameters
ProjectionProjection matrix
ResidualsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aProjectionProjection from local to measurement system (derivative of measurement vs local parameters)
[in]aResidualsMeasurement residuals
[in]aPrecisionMeasurement precision (matrix)
[in]minPrecisionMinimal precision to accept measurement

Definition at line 290 of file GblPoint.h.

Referenced by example1(), example2(), example3(), exampleDc(), and exampleSit().

◆ addMeasurement() [2/4]

template<typename Projection , typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void gbl::GblPoint::addMeasurement ( const Eigen::MatrixBase< Projection > &  aProjection,
const Eigen::MatrixBase< Residuals > &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision,
double  minPrecision = 0. 
)

Add a measurement to a point.

Add measurement (in measurement system) with diagonal precision (inverse covariance) matrix. ((up to) 2D: position, 4D: slope+position, 5D: curvature+slope+position)

Template Parameters
ProjectionProjection matrix
ResidualsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aProjectionProjection from local to measurement system (derivative of measurement vs local parameters)
[in]aResidualsMeasurement residuals
[in]aPrecisionMeasurement precision (vector with diagonal)
[in]minPrecisionMinimal precision to accept measurement

◆ addMeasurement() [3/4]

template<typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * >
void gbl::GblPoint::addMeasurement ( const Eigen::MatrixBase< Residuals > &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision,
double  minPrecision = 0. 
)

Add a measurement to a point.

Add measurement in local system with arbitrary precision (inverse covariance) matrix. Will be diagonalized. ((up to) 2D: position, 4D: slope+position, 5D: curvature+slope+position)

Template Parameters
ResidualsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aResidualsMeasurement residuals
[in]aPrecisionMeasurement precision (matrix)
[in]minPrecisionMinimal precision to accept measurement

Definition at line 332 of file GblPoint.h.

◆ addMeasurement() [4/4]

template<typename Residuals , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void gbl::GblPoint::addMeasurement ( const Eigen::MatrixBase< Residuals > &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision,
double  minPrecision = 0. 
)

Add a measurement to a point.

Add measurement in local system with diagonal precision (inverse covariance) matrix. ((up to) 2D: position, 4D: slope+position, 5D: curvature+slope+position)

Template Parameters
ResidualsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aResidualsMeasurement residuals
[in]aPrecisionMeasurement precision (vector with diagonal)
[in]minPrecisionMinimal precision to accept measurement

◆ addNextJacobian()

void gbl::GblPoint::addNextJacobian ( const Matrix5d aJac)
private

Define jacobian to next scatterer (by GBLTrajectory constructor)

Parameters
[in]aJacJacobian

Definition at line 460 of file GblPoint.cpp.

References nextJacobian.

Referenced by gbl::GblTrajectory::calcJacobians().

◆ addPrevJacobian()

void gbl::GblPoint::addPrevJacobian ( const Matrix5d aJac)
private

Define jacobian to previous scatterer (by GBLTrajectory constructor)

Parameters
[in]aJacJacobian

Definition at line 445 of file GblPoint.cpp.

References prevJacobian.

◆ addScatterer() [1/2]

template<typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * >
void gbl::GblPoint::addScatterer ( const Eigen::Vector2d &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision 
)

Add a (thin) scatterer to a point.

Add scatterer with arbitrary precision (inverse covariance) matrix. Will be diagonalized. Changes local track direction.

The precision matrix for the local slopes is defined by the angular scattering error theta_0 and the scalar products c_1, c_2 of the offset directions in the local frame with the track direction:

       (1 - c_1*c_1 - c_2*c_2)   |  1 - c_1*c_1     - c_1*c_2  |
  P =  ----------------------- * |                             |
           theta_0*theta_0       |    - c_1*c_2   1 - c_2*c_2  |
Template Parameters
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aResidualsScatterer residuals
[in]aPrecisionScatterer precision (full matrix)

Definition at line 366 of file GblPoint.h.

Referenced by example1(), example2(), example3(), exampleDc(), and exampleSit().

◆ addScatterer() [2/2]

template<typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
void gbl::GblPoint::addScatterer ( const Eigen::Vector2d &  aResiduals,
const Eigen::MatrixBase< Precision > &  aPrecision 
)

Add a (thin) scatterer to a point.

Add scatterer with diagonal precision (inverse covariance) matrix. Changes local track direction.

The precision matrix for the local slopes is defined by the angular scattering error theta_0 and the scalar products c_1, c_2 of the offset directions in the local frame with the track direction:

       (1 - c_1*c_1 - c_2*c_2)   |  1 - c_1*c_1     - c_1*c_2  |
  P =  ----------------------- * |                             |
           theta_0*theta_0       |    - c_1*c_2   1 - c_2*c_2  |
Template Parameters
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aResidualsScatterer residuals
[in]aPrecisionScatterer precision (vector with diagonal)

◆ getDerivatives()

void gbl::GblPoint::getDerivatives ( int  aDirection,
Eigen::Matrix2d &  matW,
Eigen::Matrix2d &  matWJ,
Eigen::Vector2d &  vecWd 
) const

Retrieve derivatives of local track model.

Linearized track model: F_u(q/p,u',u) = J*u + S*u' + d*q/p, W is inverse of S, negated for backward propagation.

Parameters
[in]aDirectionPropagation direction (>0 forward, else backward)
[out]matWW
[out]matWJW*J
[out]vecWdW*d
Exceptions
std::overflow_error: matrix S is singular.

Definition at line 474 of file GblPoint.cpp.

References nextJacobian, and prevJacobian.

Referenced by gbl::GblTrajectory::getFitToKinkJacobian(), and gbl::GblTrajectory::getFitToLocalJacobian().

◆ getGlobalDerivatives()

void gbl::GblPoint::getGlobalDerivatives ( Eigen::MatrixXd &  aDerivatives) const

Retrieve global derivatives from a point.

Parameters
[out]aDerivativesGlobal derivatives

Definition at line 390 of file GblPoint.cpp.

References globalDerivatives.

◆ getGlobalLabels()

void gbl::GblPoint::getGlobalLabels ( std::vector< int > &  aLabels) const

Retrieve global derivatives labels from a point.

Parameters
[out]aLabelsGlobal labels

Definition at line 382 of file GblPoint.cpp.

References globalLabels.

◆ getGlobalLabelsAndDerivatives()

void gbl::GblPoint::getGlobalLabelsAndDerivatives ( unsigned int  aRow,
std::vector< int > &  aLabels,
std::vector< double > &  aDerivatives 
) const

Retrieve global derivatives from a point for a single row.

Parameters
[in]aRowRow number
[out]aLabelsGlobal labels
[out]aDerivativesGlobal derivatives

Definition at line 400 of file GblPoint.cpp.

References globalDerivatives, and globalLabels.

◆ getLabel()

unsigned int gbl::GblPoint::getLabel ( ) const

Retrieve label of point.

Definition at line 419 of file GblPoint.cpp.

References theLabel.

◆ getLocalDerivatives()

const MatrixXd & gbl::GblPoint::getLocalDerivatives ( ) const

Retrieve local derivatives from a point.

Definition at line 340 of file GblPoint.cpp.

References localDerivatives.

◆ getMeasPrecMin()

double gbl::GblPoint::getMeasPrecMin ( ) const

get precision cutoff.

Returns
minimal measurement precision (for usage)

Definition at line 194 of file GblPoint.cpp.

References measPrecMin.

◆ getMeasTransformation()

void gbl::GblPoint::getMeasTransformation ( Eigen::MatrixXd &  aTransformation) const

Get measurement transformation (from diagonalization).

Parameters
[out]aTransformationTransformation matrix

Definition at line 216 of file GblPoint.cpp.

References measDim, measTransformation, and transFlag.

◆ getMeasurement()

void gbl::GblPoint::getMeasurement ( Matrix5d aProjection,
Vector5d aResiduals,
Vector5d aPrecision 
) const

Retrieve measurement of a point.

Parameters
[out]aProjectionProjection from (diagonalized) measurement to local system
[out]aResidualsMeasurement residuals
[out]aPrecisionMeasurement precision (diagonal)

Definition at line 204 of file GblPoint.cpp.

References measDim, measPrecision, measProjection, and measResiduals.

◆ getNumGlobals()

unsigned int gbl::GblPoint::getNumGlobals ( ) const

Retrieve number of global derivatives from a point.

Definition at line 374 of file GblPoint.cpp.

References globalDerivatives.

◆ getNumLocals()

unsigned int gbl::GblPoint::getNumLocals ( ) const

Retrieve number of local derivatives from a point.

Definition at line 335 of file GblPoint.cpp.

References localDerivatives.

◆ getOffset()

int gbl::GblPoint::getOffset ( ) const

Retrieve offset for point.

Definition at line 432 of file GblPoint.cpp.

References theOffset.

Referenced by gbl::GblTrajectory::getFitToKinkJacobian(), and gbl::GblTrajectory::getFitToLocalJacobian().

◆ getP2pJacobian()

const Matrix5d & gbl::GblPoint::getP2pJacobian ( ) const

Retrieve point-to-(previous)point jacobian.

Definition at line 437 of file GblPoint.cpp.

References p2pJacobian.

Referenced by gbl::GblTrajectory::calcJacobians().

◆ getScatterer()

void gbl::GblPoint::getScatterer ( Eigen::Matrix2d &  aTransformation,
Eigen::Vector2d &  aResiduals,
Eigen::Vector2d &  aPrecision 
) const

Retrieve scatterer of a point.

Parameters
[out]aTransformationScatterer transformation from diagonalization
[out]aResidualsScatterer residuals
[out]aPrecisionScatterer precision (diagonal)

Definition at line 290 of file GblPoint.cpp.

References scatPrecision, scatResiduals, and scatTransformation.

◆ getScatTransformation()

void gbl::GblPoint::getScatTransformation ( Eigen::Matrix2d &  aTransformation) const

Get scatterer transformation (from diagonalization).

Parameters
[out]aTransformationTransformation matrix

Definition at line 301 of file GblPoint.cpp.

References scatFlag, and scatTransformation.

◆ hasMeasurement()

unsigned int gbl::GblPoint::hasMeasurement ( ) const

Check for measurement at a point.

Get dimension of measurement (0 = none).

Returns
measurement dimension

Definition at line 186 of file GblPoint.cpp.

References measDim.

◆ hasScatterer()

bool gbl::GblPoint::hasScatterer ( ) const

Check for scatterer at a point.

Definition at line 280 of file GblPoint.cpp.

References scatFlag.

◆ operator=() [1/2]

GblPoint& gbl::GblPoint::operator= ( const GblPoint )
default

◆ operator=() [2/2]

GblPoint& gbl::GblPoint::operator= ( GblPoint &&  )
default

◆ printPoint()

void gbl::GblPoint::printPoint ( unsigned int  level = 0) const

◆ setLabel()

void gbl::GblPoint::setLabel ( unsigned int  aLabel)
private

Define label of point (by GBLTrajectory constructor)

Parameters
[in]aLabelLabel identifying point

Definition at line 414 of file GblPoint.cpp.

References theLabel.

◆ setOffset()

void gbl::GblPoint::setOffset ( int  anOffset)
private

Define offset for point (by GBLTrajectory constructor)

Parameters
[in]anOffsetOffset number

Definition at line 427 of file GblPoint.cpp.

References theOffset.

Friends And Related Function Documentation

◆ GblTrajectory

friend class GblTrajectory
friend

Definition at line 259 of file GblPoint.h.

Member Data Documentation

◆ globalDerivatives

Eigen::MatrixXd gbl::GblPoint::globalDerivatives
private

Derivatives of measurement vs additional global (MP-II) parameters.

Definition at line 285 of file GblPoint.h.

Referenced by addGlobals(), getGlobalDerivatives(), getGlobalLabelsAndDerivatives(), getNumGlobals(), and printPoint().

◆ globalLabels

std::vector<int> gbl::GblPoint::globalLabels
private

Labels of global (MP-II) derivatives.

Definition at line 284 of file GblPoint.h.

Referenced by addGlobals(), getGlobalLabels(), getGlobalLabelsAndDerivatives(), and printPoint().

◆ localDerivatives

Eigen::MatrixXd gbl::GblPoint::localDerivatives
private

Derivatives of measurement vs additional local (fit) parameters.

Definition at line 283 of file GblPoint.h.

Referenced by addLocals(), getLocalDerivatives(), getNumLocals(), and printPoint().

◆ measDim

unsigned int gbl::GblPoint::measDim
private

Dimension of measurement (1-5), 0 indicates absence of measurement.

Definition at line 270 of file GblPoint.h.

Referenced by addGlobals(), addLocals(), getMeasTransformation(), getMeasurement(), hasMeasurement(), and printPoint().

◆ measPrecision

Vector5d gbl::GblPoint::measPrecision
private

Measurement precision (diagonal of inverse covariance matrix)

Definition at line 275 of file GblPoint.h.

Referenced by getMeasurement(), and printPoint().

◆ measPrecMin

double gbl::GblPoint::measPrecMin
private

Minimal measurement precision (for usage)

Definition at line 271 of file GblPoint.h.

Referenced by getMeasPrecMin(), and printPoint().

◆ measProjection

Matrix5d gbl::GblPoint::measProjection
private

Projection from measurement to local system.

Definition at line 272 of file GblPoint.h.

Referenced by getMeasurement(), and printPoint().

◆ measResiduals

Vector5d gbl::GblPoint::measResiduals
private

Measurement residuals.

Definition at line 274 of file GblPoint.h.

Referenced by getMeasurement(), and printPoint().

◆ measTransformation

Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor , 5, 5> gbl::GblPoint::measTransformation
private

Transformation of diagonalization (of meas. precision matrix)

Definition at line 278 of file GblPoint.h.

Referenced by addGlobals(), addLocals(), and getMeasTransformation().

◆ nextJacobian

Matrix5d gbl::GblPoint::nextJacobian
private

Jacobian to next scatterer (or last measurement)

Definition at line 269 of file GblPoint.h.

Referenced by addNextJacobian(), getDerivatives(), and printPoint().

◆ p2pJacobian

Matrix5d gbl::GblPoint::p2pJacobian
private

Point-to-point jacobian from previous point.

Definition at line 267 of file GblPoint.h.

Referenced by getP2pJacobian(), and printPoint().

◆ prevJacobian

Matrix5d gbl::GblPoint::prevJacobian
private

Jacobian to previous scatterer (or first measurement)

Definition at line 268 of file GblPoint.h.

Referenced by addPrevJacobian(), getDerivatives(), and printPoint().

◆ scatFlag

bool gbl::GblPoint::scatFlag
private

Scatterer present?

Definition at line 279 of file GblPoint.h.

Referenced by getScatTransformation(), hasScatterer(), and printPoint().

◆ scatPrecision

Eigen::Vector2d gbl::GblPoint::scatPrecision
private

Scattering precision (diagonal of inverse covariance matrix)

Definition at line 282 of file GblPoint.h.

Referenced by getScatterer(), and printPoint().

◆ scatResiduals

Eigen::Vector2d gbl::GblPoint::scatResiduals
private

Scattering residuals (initial kinks if iterating)

Definition at line 281 of file GblPoint.h.

Referenced by getScatterer(), and printPoint().

◆ scatTransformation

Eigen::Matrix2d gbl::GblPoint::scatTransformation
private

Transformation of diagonalization (of scat. precision matrix)

Definition at line 280 of file GblPoint.h.

Referenced by getScatterer(), and getScatTransformation().

◆ theLabel

unsigned int gbl::GblPoint::theLabel
private

Label identifying point.

Definition at line 265 of file GblPoint.h.

Referenced by getLabel(), printPoint(), and setLabel().

◆ theOffset

int gbl::GblPoint::theOffset
private

Offset number at point if not negative (else interpolation needed)

Definition at line 266 of file GblPoint.h.

Referenced by getOffset(), printPoint(), and setOffset().

◆ transFlag

bool gbl::GblPoint::transFlag
private

Transformation exists?

Definition at line 276 of file GblPoint.h.

Referenced by addGlobals(), addLocals(), getMeasTransformation(), and printPoint().


The documentation for this class was generated from the following files: