Belle II Software  release-08-01-10
GblTrajectory Class Reference

GBL trajectory. More...

#include <GblTrajectory.h>

Collaboration diagram for GblTrajectory:

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< GblPoint > &aPointList, unsigned int aLabel, const TMatrixDSym &aSeed, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
 Create new (simple) trajectory from list of points with external seed. More...
 
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &aPointaAndTransList)
 Create new composed trajectory from list of points and transformations. More...
 
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &aPointaAndTransList, const TMatrixD &extDerivatives, const TVectorD &extMeasurements, const TVectorD &extPrecisions)
 Create new composed trajectory from list of points and transformations with (independent) external measurements. More...
 
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &aPointaAndTransList, const TMatrixD &extDerivatives, const TVectorD &extMeasurements, const TMatrixDSym &extPrecisions)
 Create new composed trajectory from list of points and transformations with (correlated) external measurements. More...
 
bool isValid () const
 Retrieve validity of trajectory.
 
unsigned int getNumPoints () const
 Retrieve number of point from trajectory.
 
unsigned int getResults (int aSignedLabel, TVectorD &localPar, TMatrixDSym &localCov) const
 Get fit results at point. More...
 
unsigned int getMeasResults (unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
 Get residuals from fit at point for measurement. More...
 
unsigned int getScatResults (unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
 Get (kink) residuals from fit at point for scatterer. More...
 
unsigned int getLabels (std::vector< unsigned int > &aLabelList)
 Get (list of) labels of points on (simple) valid trajectory. More...
 
unsigned int getLabels (std::vector< std::vector< unsigned int > > &aLabelList)
 Get (list of lists of) labels of points on (composed) valid trajectory. More...
 
unsigned int fit (double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
 Perform fit of (valid) trajectory. More...
 
void milleOut (MilleBinary &aMille)
 Write valid trajectory to Millepede-II binary file.
 
void printTrajectory (unsigned int level=0)
 Print GblTrajectory. More...
 
void printPoints (unsigned int level=0)
 Print GblPoints on trajectory. More...
 
void printData ()
 Print GblData blocks for trajectory.
 
std::vector< GblDatagetData ()
 

Private Member Functions

std::pair< std::vector< unsigned int >, TMatrixD > getJacobian (int aSignedLabel) const
 Get jacobian for transformation from fit to track parameters at point. More...
 
void getFitToLocalJacobian (std::vector< unsigned int > &anIndex, SMatrix55 &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...
 
void getFitToKinkJacobian (std::vector< unsigned int > &anIndex, SMatrix27 &aJacobian, const GblPoint &aPoint) const
 Get jacobian for transformation from (trajectory) fit to kink 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.
 
void prepare ()
 Prepare fit for simple or composed trajectory. More...
 
void buildLinearEquationSystem ()
 Build linear equation system from data (blocks).
 
void predict ()
 Calculate predictions for all points.
 
double downWeight (unsigned int aMethod)
 Down-weight all points. More...
 
void getResAndErr (unsigned int aData, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
 Get residual and errors from data block. More...
 

Private Attributes

unsigned int numAllPoints
 Number of all points on trajectory.
 
std::vector< unsigned int > numPoints
 Number of points on (sub)trajectory.
 
unsigned int numTrajectories
 Number of trajectories (in composed trajectory)
 
unsigned int numOffsets
 Number of (points with) offsets on trajectory.
 
unsigned int numInnerTrans
 Number of inner transformations to external parameters.
 
unsigned int numCurvature
 Number of curvature parameters (0 or 1) or external parameters.
 
unsigned int numParameters
 Number of fit parameters.
 
unsigned int numLocals
 Total number of (additional) local parameters.
 
unsigned int numMeasurements
 Total number of measurements.
 
unsigned int externalPoint
 Label of external point (or 0)
 
bool constructOK
 Trajectory has been successfully constructed (ready for fit/output)
 
bool fitOK
 Trajectory has been successfully fitted (results are valid)
 
std::vector< unsigned int > theDimension
 List of active dimensions (0=u1, 1=u2) in fit.
 
std::vector< std::vector< GblPoint > > thePoints
 (list of) List of points on trajectory
 
std::vector< GblDatatheData
 List of data blocks.
 
std::vector< unsigned int > measDataIndex
 mapping points to data blocks from measurements
 
std::vector< unsigned int > scatDataIndex
 mapping points to data blocks from scatterers
 
TMatrixDSym externalSeed
 Precision (inverse covariance matrix) of external seed.
 
std::vector< TMatrixD > innerTransformations
 Transformations at innermost points of.
 
TMatrixD externalDerivatives
 
TVectorD externalMeasurements
 
TVectorD externalPrecisions
 
VVector theVector
 Vector of linear equation system.
 
BorderedBandMatrix theMatrix
 (Bordered band) matrix of linear equation system
 

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 48 of file GblTrajectory.h.

Constructor & Destructor Documentation

◆ GblTrajectory() [1/5]

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

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 147 of file GblTrajectory.cc.

148  :
149  numAllPoints(aPointList.size()), numPoints(), numOffsets(0), numInnerTrans(
150  0), numCurvature(flagCurv ? 1 : 0), numParameters(0), numLocals(
151  0), numMeasurements(0), externalPoint(0), theDimension(0), thePoints(), theData(), measDataIndex(), scatDataIndex(), externalSeed(), innerTransformations(), externalDerivatives(), externalMeasurements(), externalPrecisions() {
152 
153  if (flagU1dir)
154  theDimension.push_back(0);
155  if (flagU2dir)
156  theDimension.push_back(1);
157  // simple (single) trajectory
158  thePoints.push_back(aPointList);
159  numPoints.push_back(numAllPoints);
160  construct(); // construct trajectory
161 }
unsigned int numCurvature
Number of curvature parameters (0 or 1) or external parameters.
Definition: GblTrajectory.h:92
std::vector< unsigned int > theDimension
List of active dimensions (0=u1, 1=u2) in fit.
Definition: GblTrajectory.h:99
std::vector< unsigned int > measDataIndex
mapping points to data blocks from measurements
unsigned int numParameters
Number of fit parameters.
Definition: GblTrajectory.h:93
unsigned int numLocals
Total number of (additional) local parameters.
Definition: GblTrajectory.h:94
TMatrixDSym externalSeed
Precision (inverse covariance matrix) of external seed.
void construct()
Construct trajectory from list of points.
std::vector< unsigned int > scatDataIndex
mapping points to data blocks from scatterers
unsigned int numMeasurements
Total number of measurements.
Definition: GblTrajectory.h:95
unsigned int numOffsets
Number of (points with) offsets on trajectory.
Definition: GblTrajectory.h:90
unsigned int numAllPoints
Number of all points on trajectory.
Definition: GblTrajectory.h:87
std::vector< std::vector< GblPoint > > thePoints
(list of) List of points on trajectory
unsigned int externalPoint
Label of external point (or 0)
Definition: GblTrajectory.h:96
std::vector< GblData > theData
List of data blocks.
std::vector< TMatrixD > innerTransformations
Transformations at innermost points of.
std::vector< unsigned int > numPoints
Number of points on (sub)trajectory.
Definition: GblTrajectory.h:88
unsigned int numInnerTrans
Number of inner transformations to external parameters.
Definition: GblTrajectory.h:91

◆ GblTrajectory() [2/5]

GblTrajectory ( const std::vector< GblPoint > &  aPointList,
unsigned int  aLabel,
const TMatrixDSym &  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.

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 175 of file GblTrajectory.cc.

◆ GblTrajectory() [3/5]

GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &  aPointsAndTransList)
explicit

Create new composed trajectory from list of points and transformations.

Composed of curved trajectory in space.

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

Definition at line 198 of file GblTrajectory.cc.

◆ GblTrajectory() [4/5]

GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &  aPointsAndTransList,
const TMatrixD &  extDerivatives,
const TVectorD &  extMeasurements,
const TVectorD &  extPrecisions 
)

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

Composed of curved trajectory in space.

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

Definition at line 224 of file GblTrajectory.cc.

◆ GblTrajectory() [5/5]

GblTrajectory ( const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &  aPointsAndTransList,
const TMatrixD &  extDerivatives,
const TVectorD &  extMeasurements,
const TMatrixDSym &  extPrecisions 
)

Create new composed trajectory from list of points and transformations with (correlated) external measurements.

Composed of curved trajectory in space.

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

Definition at line 254 of file GblTrajectory.cc.

Member Function Documentation

◆ construct()

void 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 302 of file GblTrajectory.cc.

◆ defineOffsets()

void 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 344 of file GblTrajectory.cc.

◆ downWeight()

double downWeight ( unsigned int  aMethod)
private

Down-weight all points.

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

Definition at line 1022 of file GblTrajectory.cc.

◆ fit()

unsigned int fit ( double &  Chi2,
int &  Ndf,
double &  lostWeight,
std::string  optionList = "" 
)

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)
Returns
Error code (non zero value indicates failure of fit)

Definition at line 1043 of file GblTrajectory.cc.

◆ getFitToKinkJacobian()

void getFitToKinkJacobian ( std::vector< unsigned int > &  anIndex,
SMatrix27 &  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 with non zero derivatives
[out]aJacobianCorresponding transformation matrix
[in]aPointPoint to use

Definition at line 583 of file GblTrajectory.cc.

◆ getFitToLocalJacobian()

void getFitToLocalJacobian ( std::vector< unsigned int > &  anIndex,
SMatrix55 &  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 with non 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 488 of file GblTrajectory.cc.

◆ getJacobian()

std::pair< std::vector< unsigned int >, TMatrixD > 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.

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 412 of file GblTrajectory.cc.

◆ getLabels() [1/2]

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

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 729 of file GblTrajectory.cc.

◆ getLabels() [2/2]

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

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 711 of file GblTrajectory.cc.

◆ getMeasResults()

unsigned int getMeasResults ( unsigned int  aLabel,
unsigned int &  numData,
TVectorD &  aResiduals,
TVectorD &  aMeasErrors,
TVectorD &  aResErrors,
TVectorD &  aDownWeights 
)

Get residuals from fit at point for measurement.

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 661 of file GblTrajectory.cc.

◆ getResAndErr()

void getResAndErr ( unsigned int  aData,
double &  aResidual,
double &  aMeasError,
double &  aResError,
double &  aDownWeight 
)
private

Get residual and errors from data block.

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
[out]aResErrorError of Residual (including correlations from track fit)
[out]aDownWeightDown-Weighting factor

Definition at line 756 of file GblTrajectory.cc.

◆ getResults()

unsigned int getResults ( int  aSignedLabel,
TVectorD &  localPar,
TMatrixDSym &  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 scatterers the track direction may change in between.

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

Definition at line 631 of file GblTrajectory.cc.

◆ getScatResults()

unsigned int getScatResults ( unsigned int  aLabel,
unsigned int &  numData,
TVectorD &  aResiduals,
TVectorD &  aMeasErrors,
TVectorD &  aResErrors,
TVectorD &  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 690 of file GblTrajectory.cc.

◆ prepare()

void prepare ( )
private

Prepare fit for simple or composed trajectory.

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

Definition at line 797 of file GblTrajectory.cc.

◆ printPoints()

void printPoints ( unsigned int  level = 0)

Print GblPoints on trajectory.

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

Definition at line 1183 of file GblTrajectory.cc.

◆ printTrajectory()

void printTrajectory ( unsigned int  level = 0)

Print GblTrajectory.

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

Definition at line 1117 of file GblTrajectory.cc.


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