8 #ifndef GBLTRAJECTORY_H_ 9 #define GBLTRAJECTORY_H_ 16 #include "TMatrixDSymEigen.h" 28 GblTrajectory(
const std::vector<GblPoint> &aPointList,
bool flagCurv =
true,
29 bool flagU1dir =
true,
bool flagU2dir =
true);
30 GblTrajectory(
const std::vector<GblPoint> &aPointList,
unsigned int aLabel,
31 const TMatrixDSym &aSeed,
bool flagCurv =
true,
bool flagU1dir =
32 true,
bool flagU2dir =
true);
34 const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointaAndTransList);
36 const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointaAndTransList,
37 const TMatrixD &extDerivatives,
const TVectorD &extMeasurements,
38 const TVectorD &extPrecisions);
42 unsigned int getResults(
int aSignedLabel, TVectorD &localPar,
43 TMatrixDSym &localCov)
const;
44 unsigned int getMeasResults(
unsigned int aLabel,
unsigned int &numRes,
45 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
46 TVectorD &aDownWeights);
47 unsigned int getScatResults(
unsigned int aLabel,
unsigned int &numRes,
48 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
49 TVectorD &aDownWeights);
50 void getLabels(std::vector<unsigned int> &aLabelList);
51 void getLabels(std::vector<std::vector< unsigned int> > &aLabelList);
52 unsigned int fit(
double &
Chi2,
int &Ndf,
double &lostWeight,
86 std::pair<std::vector<unsigned int>, TMatrixD>
getJacobian(
87 int aSignedLabel)
const;
90 unsigned int nJacobian = 1)
const;
100 void getResAndErr(
unsigned int aData,
double &aResidual,
101 double &aMeadsError,
double &aResError,
double &aDownWeight);
std::vector< unsigned int > scatDataIndex
mapping points to data blocks from scatterers
void construct()
Construct trajectory from list of points.
std::pair< std::vector< unsigned int >, TMatrixD > getJacobian(int aSignedLabel) const
Get jacobian for transformation from fit to track parameters at point.
void buildLinearEquationSystem()
Build linear equation system from data (blocks).
void printTrajectory(unsigned int level=0)
Print GblTrajectory.
TVectorD externalPrecisions
unsigned int numParameters
Number of fit parameters.
bool constructOK
Trajectory has been successfully constructed (ready for fit/output)
std::vector< unsigned int > measDataIndex
mapping points to data blocks from measurements
bool fitOK
Trajectory has been successfully fitted (results are valid)
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...
ROOT::Math::SMatrix< double, 2, 7 > SMatrix27
(Symmetric) Bordered Band Matrix.
std::vector< unsigned int > theDimension
List of active dimensions (0=u1, 1=u2) in fit.
TMatrixDSym externalSeed
Precision (inverse covariance matrix) of external seed.
void predict()
Calculate predictions for all points.
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.
ROOT::Math::SMatrix< double, 5, 5 > SMatrix55
unsigned int numOffsets
Number of (points with) offsets on trajectory.
bool isValid() const
Retrieve validity of trajectory.
void getLabels(std::vector< unsigned int > &aLabelList)
Get (list of) labels of points on (simple) trajectory.
unsigned int externalPoint
Label of external point (or 0)
unsigned int numInnerTrans
Number of inner transformations to external parameters.
void calcJacobians()
Calculate Jacobians to previous/next scatterer from point to point ones.
unsigned int numMeasurements
Total number of measurements.
void defineOffsets()
Define offsets from list of points.
void prepare()
Prepare fit for simple or composed trajectory.
unsigned int getMeasResults(unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
Get residuals at point from measurement.
void printData()
Print GblData blocks for trajectory.
Namespace for the general broken lines package.
double downWeight(unsigned int aMethod)
Down-weight all points.
std::vector< unsigned int > numPoints
Number of points on (sub)trajectory.
std::vector< TMatrixD > innerTransformations
Transformations at innermost points of.
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
Perform fit of trajectory.
void milleOut(MilleBinary &aMille)
Write valid trajectory to Millepede-II binary file.
std::vector< GblData > theData
List of data blocks.
TMatrixD externalDerivatives
std::vector< std::vector< GblPoint > > thePoints
(list of) List of points on trajectory
BorderedBandMatrix theMatrix
(Bordered band) matrix of linear equation system
VVector theVector
Vector of linear equation system.
unsigned int numTrajectories
Number of trajectories (in composed trajectory)
unsigned int getScatResults(unsigned int aLabel, unsigned int &numRes, TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors, TVectorD &aDownWeights)
Get (kink) residuals at point from scatterer.
Simple Vector based on std::vector<double>
GblTrajectory(const std::vector< GblPoint > &aPointList, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
Create new (simple) trajectory from list of points.
unsigned int numLocals
Total number of (additional) local parameters.
Millepede-II (binary) record.
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
unsigned int numAllPoints
Number of all points on trajectory.
void printPoints(unsigned int level=0)
Print GblPoints on trajectory.
void getResAndErr(unsigned int aData, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
Get residual and errors from data block.
unsigned int getResults(int aSignedLabel, TVectorD &localPar, TMatrixDSym &localCov) const
Get fit results at point.
TVectorD externalMeasurements
unsigned int numCurvature
Number of curvature parameters (0 or 1) or external parameters.