32 #ifndef GBLTRAJECTORY_H_
33 #define GBLTRAJECTORY_H_
52 GblTrajectory(
const std::vector<GblPoint> &aPointList,
bool flagCurv =
true,
53 bool flagU1dir =
true,
bool flagU2dir =
true);
55 const std::vector<std::pair<std::vector<GblPoint>, Eigen::MatrixXd> >& aPointsAndTransList);
70 template<
typename Seed>
71 GblTrajectory(
const std::vector<GblPoint> &aPointList,
unsigned int aLabel,
72 const Eigen::MatrixBase<Seed>& aSeed,
bool flagCurv =
true,
73 bool flagU1dir =
true,
bool flagU2dir =
true);
87 template<
typename Derivatives,
typename Measurements,
typename Precision,
88 typename std::enable_if<(Precision::ColsAtCompileTime != 1)>::type* =
91 const std::vector<std::pair<std::vector<GblPoint>, Eigen::MatrixXd> >& aPointsAndTransList,
92 const Eigen::MatrixBase<Derivatives>& extDerivatives,
93 const Eigen::MatrixBase<Measurements>& extMeasurements,
94 const Eigen::MatrixBase<Precision>& extPrecisions);
109 template<
typename Derivatives,
typename Measurements,
typename Precision,
110 typename std::enable_if<(Precision::ColsAtCompileTime == 1)>::type* =
113 const std::vector<std::pair<std::vector<GblPoint>, Eigen::MatrixXd> >& aPointsAndTransList,
114 const Eigen::MatrixBase<Derivatives>& extDerivatives,
115 const Eigen::MatrixBase<Measurements>& extMeasurements,
116 const Eigen::MatrixBase<Precision>& extPrecisions);
118 #ifdef GBL_EIGEN_SUPPORT_ROOT
120 GblTrajectory(
const std::vector<GblPoint> &aPointList,
unsigned int aLabel,
121 const TMatrixDSym &aSeed,
bool flagCurv =
true,
bool flagU1dir =
122 true,
bool flagU2dir =
true);
123 GblTrajectory(
const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList);
124 GblTrajectory(
const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList,
125 const TMatrixD &extDerivatives,
const TVectorD &extMeasurements,
126 const TVectorD &extPrecisions);
127 GblTrajectory(
const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList,
128 const TMatrixD &extDerivatives,
const TVectorD &extMeasurements,
129 const TMatrixDSym &extPrecisions);
134 unsigned int getResults(
int aSignedLabel, Eigen::VectorXd &localPar,
135 Eigen::MatrixXd &localCov)
const;
136 unsigned int getMeasResults(
unsigned int aLabel,
unsigned int &numData,
137 Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors,
138 Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights);
139 unsigned int getScatResults(
unsigned int aLabel,
unsigned int &numData,
140 Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors,
141 Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights);
142 #ifdef GBL_EIGEN_SUPPORT_ROOT
144 unsigned int getResults(
int aSignedLabel, TVectorD &localPar,
145 TMatrixDSym &localCov)
const;
146 unsigned int getMeasResults(
unsigned int aLabel,
unsigned int &numData,
147 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
148 TVectorD &aDownWeights);
149 unsigned int getScatResults(
unsigned int aLabel,
unsigned int &numData,
150 TVectorD &aResiduals, TVectorD &aMeasErrors, TVectorD &aResErrors,
151 TVectorD &aDownWeights);
153 unsigned int getLabels(std::vector<unsigned int> &aLabelList)
const;
155 std::vector<std::vector<unsigned int> > &aLabelList)
const;
156 unsigned int fit(
double &Chi2,
int &Ndf,
double &lostWeight,
157 const std::string& optionList =
"",
unsigned int aLabel = 0);
196 std::pair<std::vector<unsigned int>, Eigen::MatrixXd>
getJacobian(
197 int aSignedLabel)
const;
200 unsigned int nJacobian = 1)
const;
210 void getResAndErr(
unsigned int aData,
bool used,
double &aResidual,
211 double &aMeadsError,
double &aResError,
double &aDownWeight);
214 template<
typename Seed>
216 unsigned int aLabel,
const Eigen::MatrixBase<Seed>& aSeed,
217 bool flagCurv,
bool flagU1dir,
bool flagU2dir) :
218 numAllPoints(aPointList.size()), numPoints(), numOffsets(0), numInnerTransformations(
219 0), numInnerTransOffsets(0), numCurvature(flagCurv ? 1 : 0), numParameters(
220 0), numLocals(0), numMeasurements(0), externalPoint(aLabel), skippedMeasLabel(
221 0), maxNumGlobals(0), theDimension(0), thePoints(), theData(), measDataIndex(), scatDataIndex(), externalSeed(
222 aSeed), innerTransformations(), externalDerivatives(), externalMeasurements(), externalPrecisions() {
234 template<
typename Derivatives,
typename Measurements,
typename Precision,
235 typename std::enable_if<(Precision::ColsAtCompileTime != 1)>::type*>
237 const std::vector<std::pair<std::vector<GblPoint>, Eigen::MatrixXd> >& aPointsAndTransList,
238 const Eigen::MatrixBase<Derivatives>& extDerivatives,
239 const Eigen::MatrixBase<Measurements>& extMeasurements,
240 const Eigen::MatrixBase<Precision>& extPrecisions) :
241 numAllPoints(), numPoints(), numOffsets(0), numInnerTransformations(
242 aPointsAndTransList.size()), numParameters(0), numLocals(0), numMeasurements(
243 0), externalPoint(0), skippedMeasLabel(0), maxNumGlobals(0), theDimension(
244 0), thePoints(), theData(), measDataIndex(), scatDataIndex(), externalSeed(), innerTransformations() {
246 static_assert(
static_cast<int>(Measurements::ColsAtCompileTime) == 1,
"GblTrajectory: cols(Measurements) must be 1 (vector)");
247 static_assert(
static_cast<int>(Measurements::RowsAtCompileTime) ==
static_cast<int>(Derivatives::RowsAtCompileTime),
"GblTrajectory: rows(Measurements) and rows(Derivatives) must be equal");
248 static_assert(
static_cast<int>(Measurements::RowsAtCompileTime) ==
static_cast<int>(Precision::RowsAtCompileTime),
"GblTrajectory: rows(Measurements) and rows(Precision) must be equal");
249 static_assert(
static_cast<int>(Precision::RowsAtCompileTime) ==
static_cast<int>(Precision::ColsAtCompileTime),
"GblTrajectory: rows(Precision) and cols(Precision) must be equal");
251 Eigen::SelfAdjointEigenSolver<typename Precision::PlainObject> extEigen {
254 auto extTransformation = extEigen.eigenvectors().transpose();
255 externalDerivatives.resize(extDerivatives.rows(), extDerivatives.cols());
256 externalDerivatives = extTransformation * extDerivatives;
257 externalMeasurements.resize(extMeasurements.size());
258 externalMeasurements = extTransformation * extMeasurements;
259 externalPrecisions.resize(extMeasurements.size());
260 externalPrecisions = extEigen.eigenvalues();
262 for (
unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
263 thePoints.push_back(aPointsAndTransList[iTraj].first);
264 numPoints.push_back(thePoints.back().size());
265 numAllPoints += numPoints.back();
266 innerTransformations.push_back(aPointsAndTransList[iTraj].second);
268 theDimension.push_back(0);
269 theDimension.push_back(1);
271 numInnerTransOffsets = innerTransformations[0].rows() == 5 ? 2 : 1;
273 innerTransformations[0].rows() == 5 ?
274 innerTransformations[0].cols() :
275 innerTransformations[0].cols() + numInnerTransformations;
279 template<
typename Derivatives,
typename Measurements,
typename Precision,
280 typename std::enable_if<(Precision::ColsAtCompileTime == 1)>::type*>
282 const std::vector<std::pair<std::vector<GblPoint>, Eigen::MatrixXd> >& aPointsAndTransList,
283 const Eigen::MatrixBase<Derivatives>& extDerivatives,
284 const Eigen::MatrixBase<Measurements>& extMeasurements,
285 const Eigen::MatrixBase<Precision>& extPrecisions) :
286 numAllPoints(), numPoints(), numOffsets(0), numInnerTransformations(
287 aPointsAndTransList.size()), numParameters(0), numLocals(0), numMeasurements(
288 0), externalPoint(0), skippedMeasLabel(0), maxNumGlobals(0), theDimension(
289 0), thePoints(), theData(), measDataIndex(), scatDataIndex(), externalSeed(), innerTransformations() {
290 static_assert(
static_cast<int>(Measurements::ColsAtCompileTime) == 1,
"GblTrajectory: cols(Measurements) must be 1 (vector)");
291 static_assert(
static_cast<int>(Measurements::RowsAtCompileTime) ==
static_cast<int>(Derivatives::RowsAtCompileTime),
"GblTrajectory: rows(Measurements) and rows(Derivatives) must be equal");
292 static_assert(
static_cast<int>(Measurements::RowsAtCompileTime) ==
static_cast<int>(Precision::RowsAtCompileTime),
"GblTrajectory: rows(Measurements) and rows(Precision) must be equal");
298 for (
unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
299 thePoints.push_back(aPointsAndTransList[iTraj].first);
(Symmetric) Bordered Band Matrix.
std::vector< unsigned int > numPoints
Number of points on (sub)trajectory.
std::vector< unsigned int > scatDataIndex
mapping points to data blocks from scatterers
BorderedBandMatrix theMatrix
(Bordered band) matrix of linear equation system
std::vector< std::vector< GblPoint > > thePoints
(list of) List of points on trajectory
void printPoints(unsigned int level=0) const
Print GblPoints on trajectory.
unsigned int skippedMeasLabel
Label of point with measurements skipped in fit (for unbiased residuals) (or 0)
unsigned int numTrajectories
Number of trajectories (in composed trajectory)
void construct()
Construct trajectory from list of points.
unsigned int externalPoint
Label of external point (or 0)
Eigen::MatrixXd externalSeed
Precision (inverse covariance matrix) of external seed.
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.
unsigned int numAllPoints
Number of all points on trajectory.
double downWeight(unsigned int aMethod)
Down-weight all points.
VVector theVector
Vector of linear equation system.
void calcJacobians()
Calculate Jacobians to previous/next scatterer from point to point ones.
unsigned int numCurvature
Number of curvature parameters (0 or 1) or external parameters.
Eigen::VectorXd externalPrecisions
Precisions for external measurements of composed trajectory.
std::vector< std::array< unsigned int, 5 > > innerTransLab
Labels at innermost points of composed trajectory.
Eigen::VectorXd externalMeasurements
Residuals for external measurements of composed trajectory.
void printTrajectory(unsigned int level=0) const
Print GblTrajectory.
void predict()
Calculate predictions for all points.
void defineOffsets()
Define offsets from list of points.
std::vector< Eigen::MatrixXd > innerTransformations
Transformations at innermost points of composed trajectory (from common external parameters)
unsigned int numLocals
Total number of (additional) local parameters.
void getResAndErr(unsigned int aData, bool used, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
Get residual and errors from data block.
void milleOut(MilleBinary &aMille)
Write valid trajectory to Millepede-II binary file.
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.
unsigned int numInnerTransformations
Number of inner transformations to external parameters.
unsigned int maxNumGlobals
Max. number of global labels/derivatives per point.
bool isValid() const
Retrieve validity of trajectory.
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 meas...
GblTrajectory(const std::vector< GblPoint > &aPointList, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
Create new (simple) trajectory from list of points.
void prepare()
Prepare fit for simple or composed trajectory.
unsigned int getResults(int aSignedLabel, Eigen::VectorXd &localPar, Eigen::MatrixXd &localCov) const
Get fit results at point.
Eigen::MatrixXd externalDerivatives
Derivatives for external measurements of composed trajectory.
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, const std::string &optionList="", unsigned int aLabel=0)
Perform fit of (valid) trajectory.
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.
std::vector< GblData > theData
List of data blocks.
bool constructOK
Trajectory has been successfully constructed (ready for fit/output)
unsigned int numInnerTransOffsets
Number of (points with) offsets affected by inner transformations to external parameters.
std::vector< unsigned int > measDataIndex
mapping points to data blocks from measurements
void getFitToKinkJacobian(std::array< unsigned int, 7 > &anIndex, Matrix27d &aJacobian, const GblPoint &aPoint) const
Get jacobian for transformation from (trajectory) fit to kink parameters at point.
void buildLinearEquationSystem()
Build linear equation system from data (blocks).
void printData() const
Print GblData blocks for trajectory.
unsigned int getLabels(std::vector< unsigned int > &aLabelList) const
Get (list of) labels of points on (simple) valid trajectory.
unsigned int numParameters
Number of fit parameters.
unsigned int numMeasurements
Total number of measurements.
std::vector< Eigen::MatrixXd > innerTransDer
Derivatives at innermost points of composed trajectory.
std::pair< std::vector< unsigned int >, Eigen::MatrixXd > getJacobian(int aSignedLabel) const
Get jacobian for transformation from fit to track parameters at point.
bool fitOK
Trajectory has been successfully fitted (results are valid)
unsigned int numOffsets
Number of (points with) offsets on trajectory.
std::vector< unsigned int > theDimension
List of active dimensions (0=u1, 1=u2) in fit.
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
Millepede-II (binary) record.
Simple Vector based on std::vector<double>
Namespace for the general broken lines package.
Eigen::Matrix< double, 2, 7 > Matrix27d
Eigen::Matrix< double, 5, 5 > Matrix5d