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

GBL trajectory. More...

#include <GblTrajectory.h>

Public Member Functions

 GblTrajectory (const std::vector< GblPoint > &aPointList, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
 Create new (simple) trajectory from list of points. More...
 
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &aPointsAndTransList)
 Create new composed trajectory from list of points and transformations. More...
 
template<typename Seed >
 GblTrajectory (const std::vector< GblPoint > &aPointList, unsigned int aLabel, const Eigen::MatrixBase< Seed > &aSeed, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
 Create new (simple) trajectory from list of points with external seed. More...
 
template<typename Derivatives , typename Measurements , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime !=1)>::type * = nullptr>
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &aPointsAndTransList, const Eigen::MatrixBase< Derivatives > &extDerivatives, const Eigen::MatrixBase< Measurements > &extMeasurements, const Eigen::MatrixBase< Precision > &extPrecisions)
 Create new composed trajectory from list of points and transformations with arbitrary external measurements. More...
 
template<typename Derivatives , typename Measurements , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &aPointsAndTransList, const Eigen::MatrixBase< Derivatives > &extDerivatives, const Eigen::MatrixBase< Measurements > &extMeasurements, const Eigen::MatrixBase< Precision > &extPrecisions)
 Create new composed trajectory from list of points and transformations with independent external measurements. More...
 
virtual ~GblTrajectory ()
 
bool isValid () const
 Retrieve validity of trajectory. More...
 
unsigned int getNumPoints () const
 Retrieve number of point from trajectory. More...
 
unsigned int getExtResults (Eigen::VectorXd &extPar, Eigen::MatrixXd &extCov) const
 Get fit results for external parameters. More...
 
unsigned int getResults (int aSignedLabel, Eigen::VectorXd &localPar, Eigen::MatrixXd &localCov) const
 Get fit results at point. More...
 
unsigned int getMeasResults (unsigned int aLabel, unsigned int &numData, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors, Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights)
 Get residuals from fit at point for measurement "long list". More...
 
unsigned int getMeasResults (unsigned int aLabel, unsigned int &numData, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors)
 Get residuals from fit at point for measurement "short list". More...
 
unsigned int getScatResults (unsigned int aLabel, unsigned int &numData, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors, Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights)
 Get (kink) residuals from fit at point for scatterer. More...
 
unsigned int getLabels (std::vector< unsigned int > &aLabelList) const
 Get (list of) labels of points on (simple) valid trajectory. More...
 
unsigned int getLabels (std::vector< std::vector< unsigned int > > &aLabelList) const
 Get (list of lists of) labels of points on (composed) valid trajectory. More...
 
unsigned int fit (double &Chi2, int &Ndf, double &lostWeight, const std::string &optionList="", unsigned int aLabel=0)
 Perform fit of (valid) trajectory. More...
 
void milleOut (MilleBinary &aMille)
 Write valid trajectory to Millepede-II binary file. More...
 
void printTrajectory (unsigned int level=0) const
 Print GblTrajectory. More...
 
void printPoints (unsigned int level=0) const
 Print GblPoints on trajectory. More...
 
void printData () const
 Print GblData blocks for trajectory. More...
 

Private Member Functions

std::pair< std::vector< unsigned int >, Eigen::MatrixXd > getJacobian (int aSignedLabel) const
 Get jacobian for transformation from fit to track parameters at point. More...
 
void getFitToLocalJacobian (std::array< unsigned int, 5 > &anIndex, Matrix5d &aJacobian, const GblPoint &aPoint, unsigned int measDim, unsigned int nJacobian=1) const
 Get (part of) jacobian for transformation from (trajectory) fit to track parameters at point. More...
 
unsigned int getFitToKinkJacobian (std::array< unsigned int, 9 > &anIndex, Matrix49d &aJacobian, const GblPoint &aPoint) const
 Get jacobian for transformation from (trajectory) fit to kink parameters at point. More...
 
unsigned int getFitToStepJacobian (std::array< unsigned int, 9 > &anIndex, Matrix49d &aJacobian, const GblPoint &aPoint) const
 Get jacobian for transformation from (trajectory) fit to step parameters at point. More...
 
unsigned int getFitToKinkAndStepJacobian (std::array< unsigned int, 9 > &anIndex, Matrix49d &aJacobian, const GblPoint &aPoint) const
 Get jacobian for transformation from (trajectory) fit to kink and step parameters at point. More...
 
void construct ()
 Construct trajectory from list of points. More...
 
void defineOffsets ()
 Define offsets from list of points. More...
 
void calcJacobians ()
 Calculate Jacobians to previous/next scatterer from point to point ones. More...
 
void prepare ()
 Prepare fit for simple or composed trajectory. More...
 
void buildLinearEquationSystem ()
 Build linear equation system from data (blocks). More...
 
void predict ()
 Calculate predictions for all points. More...
 
double downWeight (unsigned int aMethod)
 Down-weight all points. More...
 
void getResAndErr (unsigned int aData, bool used, double &aResidual, double &aMeasError, double &aResError, double &aDownWeight)
 Get residual and errors from data block "long list". More...
 
void getResAndErr (unsigned int aData, double &aResidual, double &aMeasError)
 Get residual and errors from data block "short list". More...
 

Private Attributes

unsigned int numAllPoints
 Number of all points on trajectory. More...
 
std::vector< unsigned int > numPoints
 Number of points on (sub)trajectory. More...
 
unsigned int numTrajectories
 Number of trajectories (in composed trajectory) More...
 
unsigned int numOffsetPoints
 Number of points with offsets on trajectory. More...
 
unsigned int numOffsets
 Number of (1D or 2D) offsets on trajectory. More...
 
unsigned int numInnerTransformations
 Number of inner transformations to external parameters. More...
 
unsigned int numInnerTransOffsets
 Number of (points with) offsets affected by inner transformations to external parameters. More...
 
unsigned int numCurvature
 Number of curvature parameters (0 or 1) or external parameters. More...
 
unsigned int numParameters
 Number of fit parameters. More...
 
unsigned int numLocals
 Total number of (additional) local parameters. More...
 
unsigned int numMeasurements
 Total number of measurements. More...
 
unsigned int externalPoint
 Label of external point (or 0) More...
 
unsigned int skippedMeasLabel
 Label of point with measurements skipped in fit (for unbiased residuals) (or 0) More...
 
unsigned int maxNumGlobals
 Max. number of global labels/derivatives per point. More...
 
bool constructOK
 Trajectory has been successfully constructed (ready for fit/output) More...
 
bool fitOK
 Trajectory has been successfully fitted (results are valid) More...
 
std::vector< unsigned int > theDimension
 List of active dimensions (0=u1, 1=u2) in fit. More...
 
std::vector< std::vector< GblPoint > > thePoints
 (list of) List of points on trajectory More...
 
std::vector< GblDatatheData
 List of data blocks. More...
 
std::vector< unsigned int > measDataIndex
 mapping points to data blocks from measurements More...
 
std::vector< unsigned int > scatDataIndex
 mapping points to data blocks from scatterers More...
 
Eigen::MatrixXd externalSeed
 Precision (inverse covariance matrix) of external seed. More...
 
std::vector< Eigen::MatrixXd > innerTransformations
 Transformations at innermost points of composed trajectory (from common external parameters) More...
 
std::vector< Eigen::MatrixXd > innerTransDer
 Derivatives at innermost points of composed trajectory. More...
 
std::vector< std::array< unsigned int, 5 > > innerTransLab
 Labels at innermost points of composed trajectory. More...
 
Eigen::MatrixXd externalDerivatives
 Derivatives for external measurements of composed trajectory. More...
 
Eigen::VectorXd externalMeasurements
 Residuals for external measurements of composed trajectory. More...
 
Eigen::VectorXd externalPrecisions
 Precisions for external measurements of composed trajectory. More...
 
VVector theVector
 Vector of linear equation system. More...
 
BorderedBandMatrix theMatrix
 (Bordered band) matrix of linear equation system More...
 

Detailed Description

GBL trajectory.

List of GblPoints ordered by arc length. Can be fitted and optionally written to MP-II binary file.

Definition at line 50 of file GblTrajectory.h.

Constructor & Destructor Documentation

◆ GblTrajectory() [1/5]

gbl::GblTrajectory::GblTrajectory ( const std::vector< GblPoint > &  aPointList,
bool  flagCurv = true,
bool  flagU1dir = true,
bool  flagU2dir = true 
)

Create new (simple) trajectory from list of points.

Curved trajectory in space (default) or without curvature (q/p) or in one plane (u-direction) only.

Parameters
[in]aPointListList of points
[in]flagCurvUse q/p
[in]flagU1dirUse in u1 direction
[in]flagU2dirUse in u2 direction

Definition at line 173 of file GblTrajectory.cpp.

References construct(), numAllPoints, numPoints, theDimension, and thePoints.

◆ GblTrajectory() [2/5]

gbl::GblTrajectory::GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &  aPointsAndTransList)

Create new composed trajectory from list of points and transformations.

Composed of curved trajectories in space.

Parameters
[in]aPointsAndTransListList containing pairs with list of points and transformation (at inner (first) point)

Definition at line 196 of file GblTrajectory.cpp.

References construct(), innerTransformations, numAllPoints, numCurvature, numInnerTransformations, numInnerTransOffsets, numPoints, theDimension, and thePoints.

◆ GblTrajectory() [3/5]

template<typename Seed >
gbl::GblTrajectory::GblTrajectory ( const std::vector< GblPoint > &  aPointList,
unsigned int  aLabel,
const Eigen::MatrixBase< Seed > &  aSeed,
bool  flagCurv = true,
bool  flagU1dir = true,
bool  flagU2dir = true 
)

Create new (simple) trajectory from list of points with external seed.

Curved trajectory in space (default) or without curvature (q/p) or in one plane (u-direction) only.

Template Parameters
SeedSeed precision matrix
Parameters
[in]aPointListList of points
[in]aLabel(Signed) label of point for external seed (<0: in front, >0: after point, slope changes at scatterer!)
[in]aSeedPrecision matrix of external seed
[in]flagCurvUse q/p
[in]flagU1dirUse in u1 direction
[in]flagU2dirUse in u2 direction

Definition at line 232 of file GblTrajectory.h.

References construct(), numAllPoints, numPoints, theDimension, and thePoints.

◆ GblTrajectory() [4/5]

template<typename Derivatives , typename Measurements , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * >
gbl::GblTrajectory::GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &  aPointsAndTransList,
const Eigen::MatrixBase< Derivatives > &  extDerivatives,
const Eigen::MatrixBase< Measurements > &  extMeasurements,
const Eigen::MatrixBase< Precision > &  extPrecisions 
)

Create new composed trajectory from list of points and transformations with arbitrary external measurements.

Composed of curved trajectories in space. The precision matrix for the external measurements is specified as matrix.

Template Parameters
DerivativesExternal derivatives
MeasurementsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aPointsAndTransListList containing pairs with list of points and transformation (at inner (first) point)
[in]extDerivativesDerivatives of external measurements vs external parameters
[in]extMeasurementsExternal measurements (residuals)
[in]extPrecisionsPrecision of external measurements (matrix)

Definition at line 254 of file GblTrajectory.h.

◆ GblTrajectory() [5/5]

template<typename Derivatives , typename Measurements , typename Precision , typename std::enable_if<(Precision::ColsAtCompileTime==1)>::type * = nullptr>
gbl::GblTrajectory::GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, Eigen::MatrixXd > > &  aPointsAndTransList,
const Eigen::MatrixBase< Derivatives > &  extDerivatives,
const Eigen::MatrixBase< Measurements > &  extMeasurements,
const Eigen::MatrixBase< Precision > &  extPrecisions 
)

Create new composed trajectory from list of points and transformations with independent external measurements.

Composed of curved trajectories in space. The (diagonal) precision matrix for the external measurements is specified as vector (containing the diagonal).

Template Parameters
DerivativesExternal derivatives
MeasurementsResiduals vector
PrecisionPrecision matrix or vector (with diagonal)
Parameters
[in]aPointsAndTransListList containing pairs with list of points and transformation (at inner (first) point)
[in]extDerivativesDerivatives of external measurements vs external parameters
[in]extMeasurementsExternal measurements (residuals)
[in]extPrecisionsPrecision of external measurements (vector with diagonal)

◆ ~GblTrajectory()

gbl::GblTrajectory::~GblTrajectory ( )
virtual

Definition at line 418 of file GblTrajectory.cpp.

Member Function Documentation

◆ buildLinearEquationSystem()

void gbl::GblTrajectory::buildLinearEquationSystem ( )
private

◆ calcJacobians()

void gbl::GblTrajectory::calcJacobians ( )
private

Calculate Jacobians to previous/next scatterer from point to point ones.

Definition at line 524 of file GblTrajectory.cpp.

References gbl::GblPoint::addNextJacobian(), gbl::GblPoint::getP2pJacobian(), numTrajectories, and thePoints.

Referenced by construct().

◆ construct()

void gbl::GblTrajectory::construct ( )
private

Construct trajectory from list of points.

Trajectory is prepared for fit or output to binary file, may consists of sub-trajectories.

Definition at line 435 of file GblTrajectory.cpp.

References calcJacobians(), constructOK, defineOffsets(), fitOK, numAllPoints, numCurvature, numInnerTransformations, numInnerTransOffsets, numLocals, numMeasurements, numOffsets, numParameters, numTrajectories, prepare(), theDimension, and thePoints.

Referenced by GblTrajectory().

◆ defineOffsets()

void gbl::GblTrajectory::defineOffsets ( )
private

Define offsets from list of points.

Define offsets at points with scatterers and first and last point. All other points need interpolation from adjacent points with offsets.

Definition at line 490 of file GblTrajectory.cpp.

References numOffsetPoints, numOffsets, numTrajectories, and thePoints.

Referenced by construct().

◆ downWeight()

double gbl::GblTrajectory::downWeight ( unsigned int  aMethod)
private

Down-weight all points.

Parameters
[in]aMethodM-estimator (1: Tukey, 2:Huber, 3:Cauchy)

Definition at line 1651 of file GblTrajectory.cpp.

References theData.

Referenced by fit().

◆ fit()

unsigned int gbl::GblTrajectory::fit ( double &  Chi2,
int &  Ndf,
double &  lostWeight,
const std::string &  optionList = "",
unsigned int  aLabel = 0 
)

Perform fit of (valid) trajectory.

Optionally iterate for outlier down-weighting. Fit may fail due to singular or not positive definite matrices (internal exceptions 1-3).

Parameters
[out]Chi2Chi2 sum (corrected for down-weighting)
[out]NdfNumber of degrees of freedom
[out]lostWeightSum of weights lost due to down-weighting
[in]optionListIterations for down-weighting (One character per iteration: t,h,c (or T,H,C) for Tukey, Huber or Cauchy function)
[in]aLabelLabel of point where to skip measurements (for unbiased residuals)
Returns
Error code (non zero value indicates failure of fit)

Definition at line 1673 of file GblTrajectory.cpp.

References buildLinearEquationSystem(), constructOK, downWeight(), fitOK, gbl::InternalMeasurement, numParameters, predict(), skippedMeasLabel, gbl::BorderedBandMatrix::solveAndInvertBorderedBand(), theData, theMatrix, and theVector.

Referenced by example1(), example2(), example3(), example4(), exampleComposedGeo(), exampleComposedKin(), exampleDc(), exampleSit(), and GblTrajectory_fit().

◆ getExtResults()

unsigned int gbl::GblTrajectory::getExtResults ( Eigen::VectorXd &  extPar,
Eigen::MatrixXd &  extCov 
) const

Get fit results for external parameters.

Get corrections and covariance matrix for external parameters.

Parameters
[out]extParCorrections for external parameters
[out]extCovCovariance for external parameters
Returns
error code (non-zero if trajectory not fitted successfully)

Definition at line 922 of file GblTrajectory.cpp.

References fitOK, gbl::BorderedBandMatrix::getBlockMatrix(), innerTransformations, numLocals, theMatrix, and theVector.

◆ getFitToKinkAndStepJacobian()

unsigned int gbl::GblTrajectory::getFitToKinkAndStepJacobian ( std::array< unsigned int, 9 > &  anIndex,
Matrix49d aJacobian,
const GblPoint aPoint 
) const
private

Get jacobian for transformation from (trajectory) fit to kink and step parameters at point.

Jacobian broken lines (q/p,..,u_i-1,u_i-,u_i+,u_i+1..) to kink (du') and step (du) parameters.

Parameters
[out]anIndexList of fit parameters (zero for zero derivatives)
[out]aJacobianCorresponding transformation matrix
[in]aPointPoint to use
Returns
Number of derivatives

Definition at line 870 of file GblTrajectory.cpp.

References gbl::GblPoint::getDerivatives(), gbl::GblPoint::getOffset(), numCurvature, numLocals, and theDimension.

Referenced by prepare().

◆ getFitToKinkJacobian()

unsigned int gbl::GblTrajectory::getFitToKinkJacobian ( std::array< unsigned int, 9 > &  anIndex,
Matrix49d aJacobian,
const GblPoint aPoint 
) const
private

Get jacobian for transformation from (trajectory) fit to kink parameters at point.

Jacobian broken lines (q/p,..,u_i-1,u_i,u_i+1..) to kink (du') parameters.

Parameters
[out]anIndexList of fit parameters (zero for zero derivatives)
[out]aJacobianCorresponding transformation matrix
[in]aPointPoint to use
Returns
Number of derivatives

Definition at line 788 of file GblTrajectory.cpp.

References gbl::GblPoint::getDerivatives(), gbl::GblPoint::getOffset(), numCurvature, numLocals, and theDimension.

Referenced by prepare().

◆ getFitToLocalJacobian()

void gbl::GblTrajectory::getFitToLocalJacobian ( std::array< unsigned int, 5 > &  anIndex,
Matrix5d aJacobian,
const GblPoint aPoint,
unsigned int  measDim,
unsigned int  nJacobian = 1 
) const
private

Get (part of) jacobian for transformation from (trajectory) fit to track parameters at point.

Jacobian broken lines (q/p,..,u_i,u_i+1..) to local (q/p,u',u) parameters.

Parameters
[out]anIndexList of fit parameters (zero for zero derivatives)
[out]aJacobianCorresponding transformation matrix
[in]aPointPoint to use
[in]measDimDimension of 'measurement' (<=2: calculate only offset part, >2: complete matrix)
[in]nJacobianDirection (0: to previous offset, 1: to next offset)

Definition at line 687 of file GblTrajectory.cpp.

References gbl::GblPoint::getDerivatives(), gbl::GblPoint::getOffset(), gbl::GblPoint::getScatDim(), numCurvature, numLocals, and theDimension.

Referenced by getJacobian(), and prepare().

◆ getFitToStepJacobian()

unsigned int gbl::GblTrajectory::getFitToStepJacobian ( std::array< unsigned int, 9 > &  anIndex,
Matrix49d aJacobian,
const GblPoint aPoint 
) const
private

Get jacobian for transformation from (trajectory) fit to step parameters at point.

Jacobian broken lines (q/p,..,u_i-1,u_i-,u_i+,u_i+1..) to step (du) parameters.

Parameters
[out]anIndexList of fit parameters (zero for zero derivatives)
[out]aJacobianCorresponding transformation matrix
[in]aPointPoint to use
Returns
Number of derivatives

Definition at line 834 of file GblTrajectory.cpp.

References gbl::GblPoint::getOffset(), numCurvature, numLocals, and theDimension.

Referenced by prepare().

◆ getJacobian()

std::pair< std::vector< unsigned int >, MatrixXd > gbl::GblTrajectory::getJacobian ( int  aSignedLabel) const
private

Get jacobian for transformation from fit to track parameters at point.

Jacobian broken lines (q/p,..,u_i,u_i+1..) to track (q/p,u',u) parameters including additional local parameters (and transformation from external parameters).

Parameters
[in]aSignedLabel(Signed) label of point for external seed (<0: in front, >0: after point, slope changes at scatterer!)
Returns
List of fit parameters with non zero derivatives and corresponding transformation matrix

Definition at line 570 of file GblTrajectory.cpp.

References getFitToLocalJacobian(), innerTransDer, innerTransLab, numCurvature, numInnerTransformations, numInnerTransOffsets, numLocals, numPoints, numTrajectories, theDimension, and thePoints.

Referenced by getResults(), and prepare().

◆ getLabels() [1/2]

unsigned int gbl::GblTrajectory::getLabels ( std::vector< std::vector< unsigned int > > &  aLabelList) const

Get (list of lists of) labels of points on (composed) valid trajectory.

Parameters
[out]aLabelListList of of lists of labels
Returns
error code (non-zero if trajectory not valid (constructed successfully))

Definition at line 1216 of file GblTrajectory.cpp.

References constructOK, numTrajectories, and thePoints.

◆ getLabels() [2/2]

unsigned int gbl::GblTrajectory::getLabels ( std::vector< unsigned int > &  aLabelList) const

Get (list of) labels of points on (simple) valid trajectory.

Parameters
[out]aLabelListList of labels (aLabelList[i] = i+1)
Returns
error code (non-zero if trajectory not valid (constructed successfully))

Definition at line 1197 of file GblTrajectory.cpp.

References constructOK, and thePoints.

◆ getMeasResults() [1/2]

unsigned int gbl::GblTrajectory::getMeasResults ( unsigned int  aLabel,
unsigned int &  numData,
Eigen::VectorXd &  aResiduals,
Eigen::VectorXd &  aMeasErrors 
)

Get residuals from fit at point for measurement "short list".

Get (diagonalized) residual, error of measurement and residual and down-weighting factor for measurement at point

Parameters
[in]aLabelLabel of point on trajectory
[out]numDataNumber of data blocks from measurement at point
[out]aResidualsMeasurements-Predictions
[out]aMeasErrorsErrors of Measurements
Returns
error code (non-zero if trajectory not fitted successfully)

Definition at line 1013 of file GblTrajectory.cpp.

References fitOK, getResAndErr(), and measDataIndex.

◆ getMeasResults() [2/2]

unsigned int gbl::GblTrajectory::getMeasResults ( unsigned int  aLabel,
unsigned int &  numData,
Eigen::VectorXd &  aResiduals,
Eigen::VectorXd &  aMeasErrors,
Eigen::VectorXd &  aResErrors,
Eigen::VectorXd &  aDownWeights 
)

Get residuals from fit at point for measurement "long list".

Get (diagonalized) residual, error of measurement and residual and down-weighting factor for measurement at point

Parameters
[in]aLabelLabel of point on trajectory
[out]numDataNumber of data blocks from measurement at point
[out]aResidualsMeasurements-Predictions
[out]aMeasErrorsErrors of Measurements
[out]aResErrorsErrors of Residuals (including correlations from track fit)
[out]aDownWeightsDown-Weighting factors
Returns
error code (non-zero if trajectory not fitted successfully)

Definition at line 985 of file GblTrajectory.cpp.

References fitOK, getResAndErr(), measDataIndex, and skippedMeasLabel.

Referenced by GblTrajectory_getMeasResults().

◆ getNumPoints()

unsigned int gbl::GblTrajectory::getNumPoints ( ) const

Retrieve number of point from trajectory.

Definition at line 427 of file GblTrajectory.cpp.

References numAllPoints.

Referenced by GblTrajectory_delete(), and GblTrajectory_getNumPoints().

◆ getResAndErr() [1/2]

void gbl::GblTrajectory::getResAndErr ( unsigned int  aData,
bool  used,
double &  aResidual,
double &  aMeasError,
double &  aResError,
double &  aDownWeight 
)
private

Get residual and errors from data block "long list".

Get residual, error of measurement and residual and down-weighting factor for (single) data block

Parameters
[in]aDataLabel of data block
[in]usedFlag for usage of data block in fit
[out]aResidualMeasurement-Prediction
[out]aMeasErrorError of Measurement
[out]aResErrorError of Residual (including correlations from track fit)
[out]aDownWeightDown-Weighting factor

Definition at line 1244 of file GblTrajectory.cpp.

References gbl::BorderedBandMatrix::getBlockMatrix(), theData, and theMatrix.

Referenced by getMeasResults(), and getScatResults().

◆ getResAndErr() [2/2]

void gbl::GblTrajectory::getResAndErr ( unsigned int  aData,
double &  aResidual,
double &  aMeasError 
)
private

Get residual and errors from data block "short list".

Get residual, error of measurement and residual and down-weighting factor for (single) data block

Parameters
[in]aDataLabel of data block
[out]aResidualMeasurement-Prediction
[out]aMeasErrorError of Measurement

Definition at line 1276 of file GblTrajectory.cpp.

References theData.

◆ getResults()

unsigned int gbl::GblTrajectory::getResults ( int  aSignedLabel,
Eigen::VectorXd &  localPar,
Eigen::MatrixXd &  localCov 
) const

Get fit results at point.

Get corrections and covariance matrix for local track and additional parameters in forward or backward direction.

The point is identified by its label (1..number(points)), the sign distinguishes the backward (facing previous point) and forward 'side' (facing next point). For (thick) scatterers the track direction (and offset) may change in between.

Parameters
[in]aSignedLabel(Signed) label of point on trajectory (<0: in front, >0: after point, slope (and offset) changes at (thick) scatterer!)
[out]localParCorrections for local parameters
[out]localCovCovariance for local parameters
Returns
error code (non-zero if trajectory not fitted successfully)

Definition at line 954 of file GblTrajectory.cpp.

References fitOK, gbl::BorderedBandMatrix::getBlockMatrix(), getJacobian(), theMatrix, and theVector.

Referenced by example3(), and GblTrajectory_getResults().

◆ getScatResults()

unsigned int gbl::GblTrajectory::getScatResults ( unsigned int  aLabel,
unsigned int &  numData,
Eigen::VectorXd &  aResiduals,
Eigen::VectorXd &  aMeasErrors,
Eigen::VectorXd &  aResErrors,
Eigen::VectorXd &  aDownWeights 
)

Get (kink) residuals from fit at point for scatterer.

Get (diagonalized) residual, error of measurement and residual and down-weighting factor for scatterering kinks at point

Parameters
[in]aLabelLabel of point on trajectory
[out]numDataNumber of data blocks from scatterer at point
[out]aResiduals(kink)Measurements-(kink)Predictions
[out]aMeasErrorsErrors of (kink)Measurements
[out]aResErrorsErrors of Residuals (including correlations from track fit)
[out]aDownWeightsDown-Weighting factors
Returns
error code (non-zero if trajectory not fitted successfully)

Definition at line 1041 of file GblTrajectory.cpp.

References fitOK, getResAndErr(), and scatDataIndex.

◆ isValid()

bool gbl::GblTrajectory::isValid ( ) const

Retrieve validity of trajectory.

Definition at line 422 of file GblTrajectory.cpp.

References constructOK.

Referenced by example3(), and GblTrajectory_isValid().

◆ milleOut()

void gbl::GblTrajectory::milleOut ( MilleBinary aMille)

Write valid trajectory to Millepede-II binary file.

Trajectory state after construction (independent of fitting) is used.

Definition at line 1730 of file GblTrajectory.cpp.

References gbl::MilleBinary::addData(), constructOK, gbl::InternalMeasurement, maxNumGlobals, theData, thePoints, and gbl::MilleBinary::writeRecord().

Referenced by exampleComposedGeo(), exampleComposedKin(), exampleDc(), exampleSit(), and GblTrajectory_milleOut().

◆ predict()

void gbl::GblTrajectory::predict ( )
private

Calculate predictions for all points.

Definition at line 1640 of file GblTrajectory.cpp.

References theData, and theVector.

Referenced by fit().

◆ prepare()

void gbl::GblTrajectory::prepare ( )
private

Prepare fit for simple or composed trajectory.

Generate data (blocks) from measurements, kinks, external seed and measurements.

Exceptions
10: inner transformation matrix with invalid number of rows (valid are 5=kinematic or 2=geometric constraint)
11: inner transformation matrix with too few columns (must be >= number of rows)
12: inner transformation matrices with varying sizes
13: too many external derivatives (must be <= number of columns of inner transformation matrix)

Definition at line 1320 of file GblTrajectory.cpp.

References gbl::GblData::addDerivatives(), externalDerivatives, gbl::ExternalMeasurement, externalMeasurements, externalPoint, externalPrecisions, gbl::ExternalSeed, externalSeed, getFitToKinkAndStepJacobian(), getFitToKinkJacobian(), getFitToLocalJacobian(), getFitToStepJacobian(), getJacobian(), innerTransDer, innerTransformations, innerTransLab, gbl::InternalKink, gbl::InternalMeasurement, maxNumGlobals, measDataIndex, numAllPoints, numCurvature, numInnerTransformations, numInnerTransOffsets, numLocals, numMeasurements, numOffsets, numTrajectories, scatDataIndex, theData, theDimension, and thePoints.

Referenced by construct().

◆ printData()

void gbl::GblTrajectory::printData ( ) const

Print GblData blocks for trajectory.

Definition at line 1851 of file GblTrajectory.cpp.

References theData.

Referenced by GblTrajectory_printData().

◆ printPoints()

void gbl::GblTrajectory::printPoints ( unsigned int  level = 0) const

Print GblPoints on trajectory.

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

Definition at line 1839 of file GblTrajectory.cpp.

References numTrajectories, and thePoints.

Referenced by GblTrajectory_printPoints().

◆ printTrajectory()

void gbl::GblTrajectory::printTrajectory ( unsigned int  level = 0) const

Member Data Documentation

◆ constructOK

bool gbl::GblTrajectory::constructOK
private

Trajectory has been successfully constructed (ready for fit/output)

Definition at line 187 of file GblTrajectory.h.

Referenced by construct(), fit(), getLabels(), isValid(), milleOut(), and printTrajectory().

◆ externalDerivatives

Eigen::MatrixXd gbl::GblTrajectory::externalDerivatives
private

Derivatives for external measurements of composed trajectory.

Definition at line 199 of file GblTrajectory.h.

Referenced by prepare(), and printTrajectory().

◆ externalMeasurements

Eigen::VectorXd gbl::GblTrajectory::externalMeasurements
private

Residuals for external measurements of composed trajectory.

Definition at line 200 of file GblTrajectory.h.

Referenced by prepare(), and printTrajectory().

◆ externalPoint

unsigned int gbl::GblTrajectory::externalPoint
private

Label of external point (or 0)

Definition at line 184 of file GblTrajectory.h.

Referenced by prepare(), and printTrajectory().

◆ externalPrecisions

Eigen::VectorXd gbl::GblTrajectory::externalPrecisions
private

Precisions for external measurements of composed trajectory.

Definition at line 201 of file GblTrajectory.h.

Referenced by prepare(), and printTrajectory().

◆ externalSeed

Eigen::MatrixXd gbl::GblTrajectory::externalSeed
private

Precision (inverse covariance matrix) of external seed.

Definition at line 194 of file GblTrajectory.h.

Referenced by prepare(), and printTrajectory().

◆ fitOK

bool gbl::GblTrajectory::fitOK
private

Trajectory has been successfully fitted (results are valid)

Definition at line 188 of file GblTrajectory.h.

Referenced by construct(), fit(), getExtResults(), getMeasResults(), getResults(), getScatResults(), and printTrajectory().

◆ innerTransDer

std::vector<Eigen::MatrixXd> gbl::GblTrajectory::innerTransDer
private

Derivatives at innermost points of composed trajectory.

Definition at line 197 of file GblTrajectory.h.

Referenced by getJacobian(), and prepare().

◆ innerTransformations

std::vector<Eigen::MatrixXd> gbl::GblTrajectory::innerTransformations
private

Transformations at innermost points of composed trajectory (from common external parameters)

Definition at line 196 of file GblTrajectory.h.

Referenced by GblTrajectory(), getExtResults(), prepare(), and printTrajectory().

◆ innerTransLab

std::vector<std::array<unsigned int, 5> > gbl::GblTrajectory::innerTransLab
private

Labels at innermost points of composed trajectory.

Definition at line 198 of file GblTrajectory.h.

Referenced by getJacobian(), and prepare().

◆ maxNumGlobals

unsigned int gbl::GblTrajectory::maxNumGlobals
private

Max. number of global labels/derivatives per point.

Definition at line 186 of file GblTrajectory.h.

Referenced by milleOut(), and prepare().

◆ measDataIndex

std::vector<unsigned int> gbl::GblTrajectory::measDataIndex
private

mapping points to data blocks from measurements

Definition at line 192 of file GblTrajectory.h.

Referenced by getMeasResults(), and prepare().

◆ numAllPoints

unsigned int gbl::GblTrajectory::numAllPoints
private

Number of all points on trajectory.

Definition at line 173 of file GblTrajectory.h.

Referenced by construct(), GblTrajectory(), getNumPoints(), prepare(), and printTrajectory().

◆ numCurvature

unsigned int gbl::GblTrajectory::numCurvature
private

Number of curvature parameters (0 or 1) or external parameters.

Definition at line 180 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), construct(), GblTrajectory(), getFitToKinkAndStepJacobian(), getFitToKinkJacobian(), getFitToLocalJacobian(), getFitToStepJacobian(), getJacobian(), and prepare().

◆ numInnerTransformations

unsigned int gbl::GblTrajectory::numInnerTransformations
private

Number of inner transformations to external parameters.

Definition at line 178 of file GblTrajectory.h.

Referenced by construct(), GblTrajectory(), getJacobian(), prepare(), and printTrajectory().

◆ numInnerTransOffsets

unsigned int gbl::GblTrajectory::numInnerTransOffsets
private

Number of (points with) offsets affected by inner transformations to external parameters.

Definition at line 179 of file GblTrajectory.h.

Referenced by construct(), GblTrajectory(), getJacobian(), prepare(), and printTrajectory().

◆ numLocals

unsigned int gbl::GblTrajectory::numLocals
private

◆ numMeasurements

unsigned int gbl::GblTrajectory::numMeasurements
private

Total number of measurements.

Definition at line 183 of file GblTrajectory.h.

Referenced by construct(), prepare(), and printTrajectory().

◆ numOffsetPoints

unsigned int gbl::GblTrajectory::numOffsetPoints
private

Number of points with offsets on trajectory.

Definition at line 176 of file GblTrajectory.h.

Referenced by defineOffsets(), and printTrajectory().

◆ numOffsets

unsigned int gbl::GblTrajectory::numOffsets
private

Number of (1D or 2D) offsets on trajectory.

Definition at line 177 of file GblTrajectory.h.

Referenced by construct(), defineOffsets(), prepare(), and printTrajectory().

◆ numParameters

unsigned int gbl::GblTrajectory::numParameters
private

Number of fit parameters.

Definition at line 181 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), construct(), fit(), and printTrajectory().

◆ numPoints

std::vector<unsigned int> gbl::GblTrajectory::numPoints
private

Number of points on (sub)trajectory.

Definition at line 174 of file GblTrajectory.h.

Referenced by GblTrajectory(), and getJacobian().

◆ numTrajectories

unsigned int gbl::GblTrajectory::numTrajectories
private

Number of trajectories (in composed trajectory)

Definition at line 175 of file GblTrajectory.h.

Referenced by calcJacobians(), construct(), defineOffsets(), getJacobian(), getLabels(), prepare(), and printPoints().

◆ scatDataIndex

std::vector<unsigned int> gbl::GblTrajectory::scatDataIndex
private

mapping points to data blocks from scatterers

Definition at line 193 of file GblTrajectory.h.

Referenced by getScatResults(), and prepare().

◆ skippedMeasLabel

unsigned int gbl::GblTrajectory::skippedMeasLabel
private

Label of point with measurements skipped in fit (for unbiased residuals) (or 0)

Definition at line 185 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), fit(), and getMeasResults().

◆ theData

std::vector<GblData> gbl::GblTrajectory::theData
private

List of data blocks.

Definition at line 191 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), downWeight(), fit(), getResAndErr(), milleOut(), predict(), prepare(), and printData().

◆ theDimension

std::vector<unsigned int> gbl::GblTrajectory::theDimension
private

◆ theMatrix

BorderedBandMatrix gbl::GblTrajectory::theMatrix
private

(Bordered band) matrix of linear equation system

Definition at line 204 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), fit(), getExtResults(), getResAndErr(), getResults(), and printTrajectory().

◆ thePoints

std::vector<std::vector<GblPoint> > gbl::GblTrajectory::thePoints
private

(list of) List of points on trajectory

Definition at line 190 of file GblTrajectory.h.

Referenced by calcJacobians(), construct(), defineOffsets(), GblTrajectory(), getJacobian(), getLabels(), milleOut(), prepare(), and printPoints().

◆ theVector

VVector gbl::GblTrajectory::theVector
private

Vector of linear equation system.

Definition at line 203 of file GblTrajectory.h.

Referenced by buildLinearEquationSystem(), fit(), getExtResults(), getResults(), predict(), and printTrajectory().


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