30 #ifndef GBLTRAJECTORY_H_
31 #define GBLTRAJECTORY_H_
38 #include "TMatrixDSymEigen.h"
50 explicit GblTrajectory(
const std::vector<GblPoint> &aPointList,
bool flagCurv =
true,
51 bool flagU1dir =
true,
bool flagU2dir =
true);
52 GblTrajectory(
const std::vector<GblPoint> &aPointList,
unsigned int aLabel,
53 const TMatrixDSym &aSeed,
bool flagCurv =
true,
bool flagU1dir =
54 true,
bool flagU2dir =
true);
56 const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointaAndTransList);
58 const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointaAndTransList,
59 const TMatrixD &extDerivatives,
const TVectorD &extMeasurements,
60 const TVectorD &extPrecisions);
62 const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointaAndTransList,
63 const TMatrixD &extDerivatives,
const TVectorD &extMeasurements,
64 const TMatrixDSym &extPrecisions);
68 unsigned int getResults(
int aSignedLabel, TVectorD &localPar,
69 TMatrixDSym &localCov)
const;
70 unsigned int getMeasResults(
unsigned int aLabel,
unsigned int &numRes,
71 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
72 TVectorD &aDownWeights);
73 unsigned int getScatResults(
unsigned int aLabel,
unsigned int &numRes,
74 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
75 TVectorD &aDownWeights);
76 unsigned int getLabels(std::vector<unsigned int> &aLabelList);
77 unsigned int getLabels(std::vector<std::vector<unsigned int> > &aLabelList);
78 unsigned int fit(
double &Chi2,
int &Ndf,
double &lostWeight,
79 std::string optionList =
"");
84 std::vector<GblData> getData() {
return theData;}
107 TMatrixD externalDerivatives;
108 TVectorD externalMeasurements;
109 TVectorD externalPrecisions;
113 std::pair<std::vector<unsigned int>, TMatrixD>
getJacobian(
114 int aSignedLabel)
const;
116 SMatrix55 &aJacobian,
const GblPoint &aPoint,
unsigned int measDim,
117 unsigned int nJacobian = 1)
const;
119 SMatrix27 &aJacobian,
const GblPoint &aPoint)
const;
127 void getResAndErr(
unsigned int aData,
double &aResidual,
128 double &aMeadsError,
double &aResError,
double &aDownWeight);
BorderedBandMatrix definition.
(Symmetric) Bordered Band Matrix.
void printPoints(unsigned int level=0)
Print GblPoints on trajectory.
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.
BorderedBandMatrix theMatrix
(Bordered band) matrix of linear equation system
void prepare()
Prepare fit for simple or composed trajectory.
void getResAndErr(unsigned int aData, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
Get residual and errors from data block.
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
Perform fit of (valid) trajectory.
std::pair< std::vector< unsigned int >, TMatrixD > getJacobian(int aSignedLabel) const
Get jacobian for transformation from fit to track parameters at point.
unsigned int numCurvature
Number of curvature parameters (0 or 1) or external parameters.
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.
VVector theVector
Vector of linear equation system.
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.
std::vector< unsigned int > measDataIndex
mapping points to data blocks from measurements
unsigned int numParameters
Number of fit parameters.
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.
bool isValid() const
Retrieve validity of trajectory.
void predict()
Calculate predictions for all points.
void printTrajectory(unsigned int level=0)
Print GblTrajectory.
unsigned int numLocals
Total number of (additional) local parameters.
void milleOut(MilleBinary &aMille)
Write valid trajectory to Millepede-II binary file.
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
TMatrixDSym externalSeed
Precision (inverse covariance matrix) of external seed.
void printData()
Print GblData blocks for trajectory.
void construct()
Construct trajectory from list of points.
void buildLinearEquationSystem()
Build linear equation system from data (blocks).
void calcJacobians()
Calculate Jacobians to previous/next scatterer from point to point ones.
GblTrajectory(const std::vector< GblPoint > &aPointList, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
Create new (simple) trajectory from list of points.
std::vector< unsigned int > scatDataIndex
mapping points to data blocks from scatterers
unsigned int numMeasurements
Total number of measurements.
bool constructOK
Trajectory has been successfully constructed (ready for fit/output)
unsigned int getLabels(std::vector< unsigned int > &aLabelList)
Get (list of) labels of points on (simple) valid trajectory.
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.
unsigned int numOffsets
Number of (points with) offsets on trajectory.
unsigned int numAllPoints
Number of all points on trajectory.
std::vector< std::vector< GblPoint > > thePoints
(list of) List of points on trajectory
double downWeight(unsigned int aMethod)
Down-weight all points.
unsigned int numTrajectories
Number of trajectories (in composed trajectory)
void defineOffsets()
Define offsets from list of points.
unsigned int externalPoint
Label of external point (or 0)
std::vector< GblData > theData
List of data blocks.
std::vector< TMatrixD > innerTransformations
Transformations at innermost points of.
unsigned int getResults(int aSignedLabel, TVectorD &localPar, TMatrixDSym &localCov) const
Get fit results at point.
std::vector< unsigned int > numPoints
Number of points on (sub)trajectory.
unsigned int numInnerTrans
Number of inner transformations to external parameters.
Millepede-II (binary) record.
Simple Vector based on std::vector<double>
Namespace for the general broken lines package.