8#include <tracking/trackFindingCDC/fitting/KarimakisMethod.h>
10#include <tracking/trackFindingCDC/fitting/EigenObservationMatrix.h>
11#include <tracking/trackFindingCDC/fitting/CDCObservations2D.h>
13#include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectory2D.h>
15#include <tracking/trackFindingCDC/geometry/UncertainPerigeeCircle.h>
16#include <tracking/trackFindingCDC/geometry/PerigeeParameters.h>
17#include <tracking/trackFindingCDC/geometry/Vector2D.h>
22using namespace TrackFindingCDC;
25 : m_lineConstrained(false)
32 size_t nObservations = observations2D.
size();
34 if (not nObservations)
return;
41 double frontX = observations2D.
getX(0);
42 double frontY = observations2D.
getY(0);
45 double backX = observations2D.
getX(nObservations - 1);
46 double backY = observations2D.
getY(nObservations - 1);
66 constexpr size_t iW = 0;
67 constexpr size_t iX = 1;
68 constexpr size_t iY = 2;
69 constexpr size_t iR2 = 3;
73 const Eigen::Matrix< double, 4, 1 >& a,
74 const Eigen::Matrix< double, 4, 4 >& c,
75 bool lineConstrained =
false)
78 if (lineConstrained) {
80 q2 = c(iX, iX) - c(iY, iY);
82 q1 = c(iX, iY) * c(iR2, iR2) - c(iX, iR2) * c(iY, iR2);
83 q2 = (c(iX, iX) - c(iY, iY)) * c(iR2, iR2) - c(iX, iR2) * c(iX, iR2) + c(iY, iR2) * c(iY, iR2);
86 double phi = 0.5 * atan2(2. * q1, q2);
88 double sinphi = sin(phi);
89 double cosphi = cos(phi);
92 if (lineConstrained) {
94 I = sinphi * a(iX) - cosphi * a(iY);
97 double kappa = (sinphi * c(iX, iR2) - cosphi * c(iY, iR2)) / c(iR2, iR2);
98 double delta = -kappa * a(iR2) + sinphi * a(iX) - cosphi * a(iY);
99 curv = 2. * kappa /
sqrt(1. - 4. * delta * kappa);
100 I = 2. * delta / (1. +
sqrt(1. - 4. * delta * kappa));
105 phi += phi > 0 ? -M_PI : M_PI;
117 const Eigen::Matrix< double, 4, 4 >& c,
118 bool lineConstrained =
false)
121 const Vector2D vecPhi = -parameters.tangential();
123 const double cosphi = vecPhi.
x();
124 const double sinphi = vecPhi.
y();
126 if (lineConstrained) {
127 double chi2 = sw * (sinphi * sinphi * c(iX, iX) - 2. * sinphi * cosphi * c(iX, iY) + cosphi * cosphi * c(iY, iY));
131 const double rho = parameters.curvature();
132 const double d = parameters.impact();
134 const double u = 1 + d * rho;
135 const double kappa = 0.5 * rho / u;
137 double chi2 = sw * u * u * (sinphi * sinphi * c(iX, iX) - 2. * sinphi * cosphi * c(iX, iY) + cosphi * cosphi * c(iY,
138 iY) - kappa * kappa * c(iR2, iR2));
147 const Eigen::Matrix< double, 4, 4 >& s,
148 bool lineConstrained =
false)
152 const double curv = parameters.curvature();
153 const double I = parameters.impact();
156 const Vector2D vecPhi = -parameters.tangential();
159 const double cosphi = vecPhi.
x();
160 const double sinphi = vecPhi.
y();
162 const double ssphi = sinphi * sinphi;
163 const double scphi = sinphi * cosphi;
164 const double ccphi = cosphi * cosphi;
166 const double rho = curv;
169 const double u = 1. + rho * d;
171 using namespace NPerigeeParameterIndices;
172 if (lineConstrained) {
173 perigeePrecision(c_Curv, c_Curv) = 0.;
174 perigeePrecision(c_Curv, c_Phi0) = 0.;
175 perigeePrecision(c_Curv, c_I) = 0.;
176 perigeePrecision(c_Phi0, c_Curv) = 0.;
177 perigeePrecision(c_I, c_Curv) = 0.;
179 perigeePrecision(c_Phi0, c_Phi0) = ccphi * s(iX, iX) + 2. * scphi * s(iX, iY) + ssphi * s(iY, iY);
180 perigeePrecision(c_Phi0, c_I) = -(cosphi * s(iX) + sinphi * s(iY));
181 perigeePrecision(c_I, c_Phi0) = perigeePrecision(c_Phi0, c_I);
182 perigeePrecision(c_I, c_I) = s(iW);
185 double sa = sinphi * s(iX) - cosphi * s(iY);
186 double sb = cosphi * s(iX) + sinphi * s(iY);
187 double sc = (ssphi - ccphi) * s(iX, iY) + scphi * (s(iX, iX) - s(iY, iY));
189 double sd = sinphi * s(iX, iR2) - cosphi * s(iY, iR2);
191 double saa = ssphi * s(iX, iX) - 2. * scphi * s(iX, iY) + ccphi * s(iY, iY);
194 double se = cosphi * s(iX, iR2) + sinphi * s(iY, iR2);
195 double sbb = ccphi * s(iX, iX) + 2. * scphi * s(iX, iY) + ssphi * s(iY, iY);
197 perigeePrecision(c_Curv, c_Curv) = 0.25 * s(iR2, iR2) - d * (sd - d * (saa + 0.5 * s(iR2) - d * (sa - 0.25 * d * s(iW))));
198 perigeePrecision(c_Curv, c_Phi0) = - u * (0.5 * se - d * (sc - 0.5 * d * sb));
199 perigeePrecision(c_Phi0, c_Curv) = perigeePrecision(c_Curv, c_Phi0);
200 perigeePrecision(c_Phi0, c_Phi0) = u * u * sbb;
202 perigeePrecision(c_Curv, c_I) = rho * (-0.5 * sd + d * saa) + 0.5 * u * s(iR2) - 0.5 * d * ((2 * u + rho * d) * sa - u * d * s(iW));
203 perigeePrecision(c_I, c_Curv) = perigeePrecision(c_Curv, c_I);
204 perigeePrecision(c_Phi0, c_I) = u * (rho * sc - u * sb);
205 perigeePrecision(c_I, c_Phi0) = perigeePrecision(c_Phi0, c_I);
206 perigeePrecision(c_I, c_I) = rho * (rho * saa - 2 * u * sa) + u * u * s(iW);
208 return perigeePrecision;
217 Eigen::Matrix< double, 4, 4> sNoL = getWXYRSumMatrix(observations2D);
220 Eigen::Matrix<double, 4, 4> aNoL = sNoL / sNoL(iW);
223 Eigen::Matrix<double, 4, 1> meansNoL = aNoL.row(iW);
226 Eigen::Matrix<double, 4, 4> cNoL = aNoL - meansNoL * meansNoL.transpose();
229 size_t ndf = observations2D.
size() - 2;
237 double chi2 = calcChi2Karimaki(resultCircle, sNoL(iW), cNoL);
Class serving as a storage of observed drift circles to present to the Riemann fitter.
double getX(int iObservation) const
Getter for the x value of the observation at the given index.
double getY(int iObservation) const
Getter for the y value of the observation at the given index.
Vector2D getCentralPoint() const
Extracts the observation center that is at the index in the middle.
void passiveMoveBy(const Vector2D &origin)
Moves all observations passively such that the given vector becomes to origin of the new coordinate s...
std::size_t size() const
Returns the number of observations stored.
Particle trajectory as it is seen in xy projection represented as a circle.
double setLocalOrigin(const Vector2D &localOrigin)
Setter for the origin of the local coordinate system.
void setLocalCircle(const UncertainPerigeeCircle &localPerigeeCircle)
Setter for the generalized circle that describes the trajectory.
void clear()
Clears all information from this trajectory.
KarimakisMethod()
Constructor setting the default constraints.
bool isLineConstrained() const
Getter for the indicator that lines should be fitted by this fitter.
UncertainPerigeeCircle fitInternal(CDCObservations2D &observations2D) const
Internal method doing the heavy work.
void update(CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
Executes the fit and updates the trajectory parameters. This may render the information in the observ...
Extension of the generalized circle also caching the perigee coordinates.
double arcLengthBetween(const Vector2D &from, const Vector2D &to) const
Calculates the arc length between two points of closest approach on the circle.
A matrix implementation to be used as an interface typ through out the track finder.
Adds an uncertainty matrix to the circle in perigee parameterisation.
void reverse()
Flips the orientation of the circle in place.
void setPerigeeCovariance(const PerigeeCovariance &perigeeCovariance)
Setter for the whole covariance matrix of the perigee parameters.
void setNDF(std::size_t ndf)
Setter for the number of degrees of freediom used in the circle fit.
void setChi2(const double chi2)
Setter for the chi square value of the circle fit.
A two dimensional vector which is equipped with functions for correct handling of orientation relate...
double x() const
Getter for the x coordinate.
double y() const
Getter for the y coordinate.
double sqrt(double a)
sqrt for double
Abstract base class for different kinds of events.
static CovarianceMatrix covarianceFromPrecision(const PrecisionMatrix &prec)
Convert the precision matrix to the corresponding covariance matrix.