GeneralBrokenLines V03-01-01
using EIGEN
Public Member Functions | Private Attributes | List of all members
gbl::GblMeasurement Class Reference

Measurement at point. More...

#include <GblMeasurement.h>

Public Member Functions

 GblMeasurement (const Eigen::MatrixXd &aProjection, const Eigen::VectorXd &aResiduals, const Eigen::MatrixXd &aPrecision, double minPrecision=0.)
 Create a measurement. More...
 
 GblMeasurement (const Eigen::VectorXd &aResiduals, const Eigen::MatrixXd &aPrecision, double minPrecision=0.)
 Create a measurement. More...
 
 GblMeasurement (const GblMeasurement &)=default
 
GblMeasurementoperator= (const GblMeasurement &)=default
 
 GblMeasurement (GblMeasurement &&)=default
 
GblMeasurementoperator= (GblMeasurement &&)=default
 
virtual ~GblMeasurement ()
 
void addLocals (const Eigen::MatrixXd &aDerivatives)
 Add local derivatives. More...
 
void addGlobals (const std::vector< int > &aLabels, const Eigen::MatrixXd &aDerivatives)
 Add global derivatives. More...
 
void setEnabled (bool flag)
 Set enabled flag. More...
 
bool isEnabled () const
 Get enabled flag. More...
 
unsigned int getMeasDim () const
 Get measurement dimension. 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...
 
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 printMeasurement (unsigned int level=0) const
 Print GblMeasurement. 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...
 

Private Attributes

bool enabled
 Enabled flag (to be used) 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...
 
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...
 

Detailed Description

Measurement at point.

User supplied measurement at point on (initial) trajectory.

Must have measurement (1D - 5D). May have:

  1. Additional local parameters (with derivatives). Fitted together with track parameters.
  2. Additional global parameters (with labels and derivatives). Not fitted, only passed on to (binary) file for fitting with Millepede-II.

Definition at line 67 of file GblMeasurement.h.

Constructor & Destructor Documentation

◆ GblMeasurement() [1/4]

gbl::GblMeasurement::GblMeasurement ( const Eigen::MatrixXd &  aProjection,
const Eigen::VectorXd &  aResiduals,
const Eigen::MatrixXd &  aPrecision,
double  minPrecision = 0. 
)

Create a measurement.

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

Parameters
[in]aProjectionProjection from local to measurement system
[in]aResidualsMeasurement residuals
[in]aPrecisionMeasurement precision (vector (with diagonal) or (full) matrix)
[in]minPrecisionMinimal precision to accept measurement

Definition at line 46 of file GblMeasurement.cpp.

References enabled, measDim, measPrecision, measPrecMin, measProjection, measResiduals, measTransformation, and transFlag.

◆ GblMeasurement() [2/4]

gbl::GblMeasurement::GblMeasurement ( const Eigen::VectorXd &  aResiduals,
const Eigen::MatrixXd &  aPrecision,
double  minPrecision = 0. 
)

Create a measurement.

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

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

Definition at line 80 of file GblMeasurement.cpp.

References enabled, measDim, measPrecision, measPrecMin, measProjection, measResiduals, measTransformation, and transFlag.

◆ GblMeasurement() [3/4]

gbl::GblMeasurement::GblMeasurement ( const GblMeasurement )
default

◆ GblMeasurement() [4/4]

gbl::GblMeasurement::GblMeasurement ( GblMeasurement &&  )
default

◆ ~GblMeasurement()

gbl::GblMeasurement::~GblMeasurement ( )
virtual

Definition at line 103 of file GblMeasurement.cpp.

Member Function Documentation

◆ addGlobals()

void gbl::GblMeasurement::addGlobals ( const std::vector< int > &  aLabels,
const Eigen::MatrixXd &  aDerivatives 
)

Add global derivatives.

Add global derivatives to measurement.

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

Definition at line 247 of file GblMeasurement.cpp.

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

◆ addLocals()

void gbl::GblMeasurement::addLocals ( const Eigen::MatrixXd &  aDerivatives)

Add local derivatives.

Add local derivatives to measurement.

Parameters
[in]aDerivativesLocal derivatives (matrix)

Definition at line 230 of file GblMeasurement.cpp.

References localDerivatives, measDim, measTransformation, and transFlag.

◆ getGlobalDerivatives()

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

Retrieve global derivatives from a point.

Parameters
[out]aDerivativesGlobal derivatives

Definition at line 402 of file GblMeasurement.cpp.

References globalDerivatives.

◆ getGlobalLabels()

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

Retrieve global derivatives labels from a point.

Parameters
[out]aLabelsGlobal labels

Definition at line 394 of file GblMeasurement.cpp.

References globalLabels.

◆ getGlobalLabelsAndDerivatives()

void gbl::GblMeasurement::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 412 of file GblMeasurement.cpp.

References globalDerivatives, and globalLabels.

◆ getLocalDerivatives()

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

Retrieve local derivatives from a point.

Definition at line 381 of file GblMeasurement.cpp.

References localDerivatives.

◆ getMeasDim()

unsigned int gbl::GblMeasurement::getMeasDim ( ) const

Get measurement dimension.

Get dimension of measurement (0 = none).

Returns
measurement dimension

Definition at line 336 of file GblMeasurement.cpp.

References measDim.

◆ getMeasPrecMin()

double gbl::GblMeasurement::getMeasPrecMin ( ) const

get precision cutoff.

Returns
minimal measurement precision (for usage)

Definition at line 344 of file GblMeasurement.cpp.

References measPrecMin.

◆ getMeasTransformation()

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

Get measurement transformation (from diagonalization).

Parameters
[out]aTransformationTransformation matrix

Definition at line 366 of file GblMeasurement.cpp.

References measDim, measTransformation, and transFlag.

◆ getMeasurement()

void gbl::GblMeasurement::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 354 of file GblMeasurement.cpp.

References measDim, measPrecision, measProjection, and measResiduals.

◆ getNumGlobals()

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

Retrieve number of global derivatives from a point.

Definition at line 386 of file GblMeasurement.cpp.

References globalDerivatives.

◆ getNumLocals()

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

Retrieve number of local derivatives from a point.

Definition at line 376 of file GblMeasurement.cpp.

References localDerivatives.

◆ isEnabled()

bool gbl::GblMeasurement::isEnabled ( ) const

Get enabled flag.

Returns
enabled flag

Definition at line 327 of file GblMeasurement.cpp.

References enabled.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ printMeasurement()

void gbl::GblMeasurement::printMeasurement ( unsigned int  level = 0) const

Print GblMeasurement.

Parameters
[in]levelprint level (0: minimum, >0: more)

Definition at line 426 of file GblMeasurement.cpp.

References enabled, globalDerivatives, globalLabels, localDerivatives, measDim, measPrecision, measPrecMin, measProjection, measResiduals, and transFlag.

◆ setEnabled()

void gbl::GblMeasurement::setEnabled ( bool  flag)

Set enabled flag.

Allows to disable measurements at points (in input list of points for trajectory constructor) to construct new trajectory with reduced number of measurements (e.g. after solution of ambiguities).

Parameters
[in]flagenabled flag

Definition at line 319 of file GblMeasurement.cpp.

References enabled.

Member Data Documentation

◆ enabled

bool gbl::GblMeasurement::enabled
private

Enabled flag (to be used)

Definition at line 115 of file GblMeasurement.h.

Referenced by GblMeasurement(), isEnabled(), printMeasurement(), and setEnabled().

◆ globalDerivatives

Eigen::MatrixXd gbl::GblMeasurement::globalDerivatives
private

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

Definition at line 127 of file GblMeasurement.h.

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

◆ globalLabels

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

Labels of global (MP-II) derivatives.

Definition at line 126 of file GblMeasurement.h.

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

◆ localDerivatives

Eigen::MatrixXd gbl::GblMeasurement::localDerivatives
private

Derivatives of measurement vs additional local (fit) parameters.

Definition at line 125 of file GblMeasurement.h.

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

◆ measDim

unsigned int gbl::GblMeasurement::measDim
private

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

Definition at line 116 of file GblMeasurement.h.

Referenced by addGlobals(), addLocals(), GblMeasurement(), getMeasDim(), getMeasTransformation(), getMeasurement(), and printMeasurement().

◆ measPrecision

Vector5d gbl::GblMeasurement::measPrecision
private

Measurement precision (diagonal of inverse covariance matrix)

Definition at line 121 of file GblMeasurement.h.

Referenced by GblMeasurement(), getMeasurement(), and printMeasurement().

◆ measPrecMin

double gbl::GblMeasurement::measPrecMin
private

Minimal measurement precision (for usage)

Definition at line 117 of file GblMeasurement.h.

Referenced by GblMeasurement(), getMeasPrecMin(), and printMeasurement().

◆ measProjection

Matrix5d gbl::GblMeasurement::measProjection
private

Projection from measurement to local system.

Definition at line 118 of file GblMeasurement.h.

Referenced by GblMeasurement(), getMeasurement(), and printMeasurement().

◆ measResiduals

Vector5d gbl::GblMeasurement::measResiduals
private

Measurement residuals.

Definition at line 120 of file GblMeasurement.h.

Referenced by GblMeasurement(), getMeasurement(), and printMeasurement().

◆ measTransformation

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

Transformation of diagonalization (of meas. precision matrix)

Definition at line 124 of file GblMeasurement.h.

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

◆ transFlag

bool gbl::GblMeasurement::transFlag
private

Transformation exists?

Definition at line 122 of file GblMeasurement.h.

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


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