GeneralBrokenLines V03-01-01
using EIGEN
GblMeasurement.h
Go to the documentation of this file.
1/*
2 * GblMeasurement.h
3 *
4 * Created on: 31 Mar 2021
5 * Author: kleinwrt
6 */
7
31#ifndef GBLMEASUREMENT_H_
32#define GBLMEASUREMENT_H_
33
34#include<iostream>
35#include<vector>
36#include <array>
37#include<math.h>
38#include "VMatrix.h"
39
40#ifdef GBL_EIGEN_SUPPORT_ROOT
41#include "TVectorD.h"
42#include "TMatrixD.h"
43#include "TMatrixDSym.h"
44#include "TMatrixDSymEigen.h"
45#endif
46
47#include "Eigen/Dense"
48
49namespace gbl {
50
51typedef Eigen::Matrix<double, 5, 1> Vector5d;
52typedef Eigen::Matrix<double, 2, 3> Matrix23d;
53typedef Eigen::Matrix<double, 2, 5> Matrix25d;
54typedef Eigen::Matrix<double, 3, 2> Matrix32d;
55typedef Eigen::Matrix<double, 5, 5> Matrix5d;
56
58
68public:
69 GblMeasurement(const Eigen::MatrixXd &aProjection,
70 const Eigen::VectorXd &aResiduals,
71 const Eigen::MatrixXd &aPrecision, double minPrecision = 0.);
72 GblMeasurement(const Eigen::VectorXd &aResiduals,
73 const Eigen::MatrixXd &aPrecision, double minPrecision = 0.);
74 GblMeasurement(const GblMeasurement&) = default;
78 virtual ~GblMeasurement();
79#ifdef GBL_EIGEN_SUPPORT_ROOT
80 // input via ROOT
81 GblMeasurement(const TMatrixD &aProjection, const TVectorD &aResiduals,
82 const TVectorD &aPrecision, double minPrecision = 0.);
83 GblMeasurement(const TMatrixD &aProjection, const TVectorD &aResiduals,
84 const TMatrixDSym &aPrecision, double minPrecision = 0.);
85 GblMeasurement(const TVectorD &aResiduals, const TVectorD &aPrecision,
86 double minPrecision = 0.);
87 GblMeasurement(const TVectorD &aResiduals, const TMatrixDSym &aPrecision,
88 double minPrecision = 0.);
89 void addLocals(const TMatrixD &aDerivatives);
90 void addGlobals(const std::vector<int> &aLabels,
91 const TMatrixD &aDerivatives);
92#endif
93
94 // input via Eigen
95 void addLocals(const Eigen::MatrixXd &aDerivatives);
96 void addGlobals(const std::vector<int> &aLabels,
97 const Eigen::MatrixXd &aDerivatives);
98 void setEnabled(bool flag);
99 bool isEnabled() const;
100 unsigned int getMeasDim() const;
101 double getMeasPrecMin() const;
102 void getMeasurement(Matrix5d &aProjection, Vector5d &aResiduals,
103 Vector5d &aPrecision) const;
104 void getMeasTransformation(Eigen::MatrixXd &aTransformation) const;
105 unsigned int getNumLocals() const;
106 const Eigen::MatrixXd& getLocalDerivatives() const;
107 unsigned int getNumGlobals() const;
108 void printMeasurement(unsigned int level = 0) const;
109 void getGlobalLabels(std::vector<int> &aLabels) const;
110 void getGlobalDerivatives(Eigen::MatrixXd &aDerivatives) const;
111 void getGlobalLabelsAndDerivatives(unsigned int aRow,
112 std::vector<int> &aLabels, std::vector<double> &aDerivatives) const;
113
114private:
115 bool enabled;
116 unsigned int measDim;
117 double measPrecMin;
119
123 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic,
124 Eigen::ColMajor /* default */, 5, 5> measTransformation;
125 Eigen::MatrixXd localDerivatives;
126 std::vector<int> globalLabels;
127 Eigen::MatrixXd globalDerivatives;
128};
129}
130
131#endif /* GBLMEASUREMENT_H_ */
VMatrix definition.
Measurement at point.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor, 5, 5 > measTransformation
Transformation of diagonalization (of meas. precision matrix)
double getMeasPrecMin() const
get precision cutoff.
unsigned int getNumGlobals() const
Retrieve number of global derivatives from a point.
Vector5d measPrecision
Measurement precision (diagonal of inverse covariance matrix)
void addLocals(const Eigen::MatrixXd &aDerivatives)
Add local derivatives.
void getGlobalLabelsAndDerivatives(unsigned int aRow, std::vector< int > &aLabels, std::vector< double > &aDerivatives) const
Retrieve global derivatives from a point for a single row.
GblMeasurement(const Eigen::MatrixXd &aProjection, const Eigen::VectorXd &aResiduals, const Eigen::MatrixXd &aPrecision, double minPrecision=0.)
Create a measurement.
void getGlobalDerivatives(Eigen::MatrixXd &aDerivatives) const
Retrieve global derivatives from a point.
Eigen::MatrixXd globalDerivatives
Derivatives of measurement vs additional global (MP-II) parameters.
bool transFlag
Transformation exists?
double measPrecMin
Minimal measurement precision (for usage)
bool enabled
Enabled flag (to be used)
GblMeasurement(const GblMeasurement &)=default
GblMeasurement & operator=(const GblMeasurement &)=default
void getMeasTransformation(Eigen::MatrixXd &aTransformation) const
Get measurement transformation (from diagonalization).
void getGlobalLabels(std::vector< int > &aLabels) const
Retrieve global derivatives labels from a point.
Matrix5d measProjection
Projection from measurement to local system.
bool isEnabled() const
Get enabled flag.
void getMeasurement(Matrix5d &aProjection, Vector5d &aResiduals, Vector5d &aPrecision) const
Retrieve measurement of a point.
void printMeasurement(unsigned int level=0) const
Print GblMeasurement.
GblMeasurement(GblMeasurement &&)=default
Eigen::MatrixXd localDerivatives
Derivatives of measurement vs additional local (fit) parameters.
unsigned int getNumLocals() const
Retrieve number of local derivatives from a point.
Vector5d measResiduals
Measurement residuals.
GblMeasurement & operator=(GblMeasurement &&)=default
std::vector< int > globalLabels
Labels of global (MP-II) derivatives.
unsigned int getMeasDim() const
Get measurement dimension.
unsigned int measDim
Dimension of measurement (1-5), 0 indicates absence of measurement.
void setEnabled(bool flag)
Set enabled flag.
const Eigen::MatrixXd & getLocalDerivatives() const
Retrieve local derivatives from a point.
void addGlobals(const std::vector< int > &aLabels, const Eigen::MatrixXd &aDerivatives)
Add global derivatives.
Namespace for the general broken lines package.
Eigen::Matrix< double, 3, 2 > Matrix32d
Eigen::Matrix< double, 2, 5 > Matrix25d
Eigen::Matrix< double, 2, 3 > Matrix23d
Eigen::Matrix< double, 5, 1 > Vector5d
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: GblData.h:46