GeneralBrokenLines  trunk_rev46
 All Classes Files Functions Variables Typedefs Pages
Public Member Functions | Private Member Functions | Private Attributes
GblTrajectory Class Reference

GBL trajectory. More...

#include <GblTrajectory.h>

List of all members.

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.
 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.
 GblTrajectory (const std::vector< std::pair< std::vector< GblPoint >, TMatrixD > > &aPointaAndTransList)
 Create new composed trajectory from list of points and transformations.
 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.
virtual ~GblTrajectory ()
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.
unsigned int getMeasResults (unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
 Get residuals at point from measurement.
unsigned int getScatResults (unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
 Get (kink) residuals at point from scatterer.
unsigned int fit (double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
 Perform fit of trajectory.
void milleOut (MilleBinary &aMille)
 Write trajectory to Millepede-II binary file.

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.
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.
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.
void construct ()
 Construct trajectory from list of points.
void defineOffsets ()
 Define offsets from list of points.
void calcJacobians ()
 Calculate Jacobians to previous/next scatterer from point to point ones.
void prepare ()
 Prepare fit for simple or composed trajectory.
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.
void getResAndErr (unsigned int aData, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
 Get residual and errors from data block.

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 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
std::vector< unsigned int > externalIndex
 List of fit parameters used by external seed.
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 23 of file GblTrajectory.h.


Constructor & Destructor Documentation

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 102 of file GblTrajectory.cpp.

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

GblTrajectory::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 130 of file GblTrajectory.cpp.

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

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

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 153 of file GblTrajectory.cpp.

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

GblTrajectory::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 179 of file GblTrajectory.cpp.

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

GblTrajectory::~GblTrajectory ( )
virtual

Definition at line 201 of file GblTrajectory.cpp.


Member Function Documentation

void GblTrajectory::buildLinearEquationSystem ( )
private

Build linear equation system from data (blocks).

Definition at line 632 of file GblTrajectory.cpp.

References BorderedBandMatrix::addBlockMatrix(), numCurvature, numLocals, numParameters, VVector::resize(), BorderedBandMatrix::resize(), theData, theMatrix, and theVector.

Referenced by fit().

void GblTrajectory::calcJacobians ( )
private

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

Definition at line 264 of file GblTrajectory.cpp.

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

Referenced by construct().

void 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 213 of file GblTrajectory.cpp.

References calcJacobians(), defineOffsets(), fitOK, numCurvature, numInnerTrans, numLocals, numMeasurements, numOffsets, numParameters, numTrajectories, prepare(), theDimension, and thePoints.

Referenced by GblTrajectory().

void 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 242 of file GblTrajectory.cpp.

References numOffsets, numTrajectories, and thePoints.

Referenced by construct().

double GblTrajectory::downWeight ( unsigned int  aMethod)
private

Down-weight all points.

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

Definition at line 873 of file GblTrajectory.cpp.

References theData.

Referenced by fit().

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

Perform fit of trajectory.

Optionally iterate for outlier down-weighting.

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 892 of file GblTrajectory.cpp.

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

Referenced by example1().

void GblTrajectory::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 486 of file GblTrajectory.cpp.

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

Referenced by prepare().

void GblTrajectory::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 392 of file GblTrajectory.cpp.

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

Referenced by getJacobian(), and prepare().

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

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 316 of file GblTrajectory.cpp.

References getFitToLocalJacobian(), numCurvature, numLocals, numPoints, numTrajectories, theDimension, and thePoints.

Referenced by getResults(), and prepare().

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

Get residuals at point from measurement.

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

Parameters:
[in]aLabelLabel of point
[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 558 of file GblTrajectory.cpp.

References fitOK, getResAndErr(), and measDataIndex.

unsigned int GblTrajectory::getNumPoints ( ) const

Retrieve number of point from trajectory.

Definition at line 205 of file GblTrajectory.cpp.

References numAllPoints.

void GblTrajectory::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 612 of file GblTrajectory.cpp.

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

Referenced by getMeasResults(), and getScatResults().

unsigned int GblTrajectory::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.

Parameters:
[in]aSignedLabel(Signed) label of point (<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 529 of file GblTrajectory.cpp.

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

unsigned int GblTrajectory::getScatResults ( unsigned int  aLabel,
unsigned int &  numData,
TVectorD &  aResiduals,
TVectorD &  aMeasErrors,
TVectorD &  aResErrors,
TVectorD &  aDownWeights 
)

Get (kink) residuals at point from scatterer.

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

Parameters:
[in]aLabelLabel of point
[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 586 of file GblTrajectory.cpp.

References fitOK, getResAndErr(), and scatDataIndex.

void GblTrajectory::milleOut ( MilleBinary aMille)

Write trajectory to Millepede-II binary file.

Definition at line 937 of file GblTrajectory.cpp.

References MilleBinary::addData(), theData, and MilleBinary::writeRecord().

void GblTrajectory::predict ( )
private

Calculate predictions for all points.

Definition at line 862 of file GblTrajectory.cpp.

References theData, and theVector.

Referenced by fit().

void GblTrajectory::prepare ( )
private

Member Data Documentation

TMatrixD GblTrajectory::externalDerivatives
private

Definition at line 71 of file GblTrajectory.h.

Referenced by prepare().

std::vector<unsigned int> GblTrajectory::externalIndex
private

List of fit parameters used by external seed.

Definition at line 67 of file GblTrajectory.h.

Referenced by prepare().

TVectorD GblTrajectory::externalMeasurements
private

Definition at line 72 of file GblTrajectory.h.

Referenced by prepare().

unsigned int GblTrajectory::externalPoint
private

Label of external point (or 0)

Definition at line 60 of file GblTrajectory.h.

Referenced by prepare().

TVectorD GblTrajectory::externalPrecisions
private

Definition at line 73 of file GblTrajectory.h.

Referenced by prepare().

TMatrixDSym GblTrajectory::externalSeed
private

Precision (inverse covariance matrix) of external seed.

Definition at line 68 of file GblTrajectory.h.

Referenced by prepare().

bool GblTrajectory::fitOK
private

Trajectory has been successfully fitted (results are valid)

Definition at line 61 of file GblTrajectory.h.

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

std::vector<TMatrixD> GblTrajectory::innerTransformations
private

Transformations at innermost points of.

Definition at line 69 of file GblTrajectory.h.

Referenced by GblTrajectory(), and prepare().

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

mapping points to data blocks from measurements

Definition at line 65 of file GblTrajectory.h.

Referenced by getMeasResults(), and prepare().

unsigned int GblTrajectory::numAllPoints
private

Number of all points on trajectory.

Definition at line 51 of file GblTrajectory.h.

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

unsigned int GblTrajectory::numCurvature
private

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

Definition at line 56 of file GblTrajectory.h.

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

unsigned int GblTrajectory::numInnerTrans
private

Number of inner transformations to external parameters.

Definition at line 55 of file GblTrajectory.h.

Referenced by construct(), and prepare().

unsigned int GblTrajectory::numLocals
private

Total number of (additional) local parameters.

Definition at line 58 of file GblTrajectory.h.

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

unsigned int GblTrajectory::numMeasurements
private

Total number of measurements.

Definition at line 59 of file GblTrajectory.h.

Referenced by construct(), and prepare().

unsigned int GblTrajectory::numOffsets
private

Number of (points with) offsets on trajectory.

Definition at line 54 of file GblTrajectory.h.

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

unsigned int GblTrajectory::numParameters
private

Number of fit parameters.

Definition at line 57 of file GblTrajectory.h.

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

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

Number of points on (sub)trajectory.

Definition at line 52 of file GblTrajectory.h.

Referenced by GblTrajectory(), and getJacobian().

unsigned int GblTrajectory::numTrajectories
private

Number of trajectories (in composed trajectory)

Definition at line 53 of file GblTrajectory.h.

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

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

mapping points to data blocks from scatterers

Definition at line 66 of file GblTrajectory.h.

Referenced by getScatResults(), and prepare().

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

List of data blocks.

Definition at line 64 of file GblTrajectory.h.

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

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

List of active dimensions (0=u1, 1=u2) in fit.

Definition at line 62 of file GblTrajectory.h.

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

BorderedBandMatrix GblTrajectory::theMatrix
private

(Bordered band) matrix of linear equation system

Definition at line 75 of file GblTrajectory.h.

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

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

(list of) List of points on trajectory

Definition at line 63 of file GblTrajectory.h.

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

VVector GblTrajectory::theVector
private

Vector of linear equation system.

Definition at line 74 of file GblTrajectory.h.

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


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