Belle II Software development
KarimakisMethod.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#include <tracking/trackFindingCDC/fitting/KarimakisMethod.h>
9
10#include <tracking/trackFindingCDC/fitting/EigenObservationMatrix.h>
11#include <tracking/trackFindingCDC/fitting/CDCObservations2D.h>
12
13#include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectory2D.h>
14
15#include <tracking/trackFindingCDC/geometry/UncertainPerigeeCircle.h>
16#include <tracking/trackFindingCDC/geometry/PerigeeParameters.h>
17#include <tracking/trackFindingCDC/geometry/Vector2D.h>
18
19#include <Eigen/Core>
20
21using namespace Belle2;
22using namespace TrackFindingCDC;
23
25 : m_lineConstrained(false)
26{
27}
28
30 CDCObservations2D& observations2D) const
31{
32 size_t nObservations = observations2D.size();
33 trajectory2D.clear();
34 if (not nObservations) return;
35
36 Vector2D ref = observations2D.getCentralPoint();
37 observations2D.passiveMoveBy(ref);
38
39 UncertainPerigeeCircle perigeeCircle = fitInternal(observations2D);
40
41 double frontX = observations2D.getX(0);
42 double frontY = observations2D.getY(0);
43 Vector2D frontPos(frontX, frontY);
44
45 double backX = observations2D.getX(nObservations - 1);
46 double backY = observations2D.getY(nObservations - 1);
47 Vector2D backPos(backX, backY);
48
49 Vector2D overPos(0, 0);
50 double totalPerps = (perigeeCircle->arcLengthBetween(frontPos, overPos) +
51 perigeeCircle->arcLengthBetween(overPos, backPos));
52
53 if (totalPerps < 0) {
54 perigeeCircle.reverse();
55 }
56
57 trajectory2D.setLocalOrigin(ref);
58 trajectory2D.setLocalCircle(perigeeCircle);
59}
60
61
62
63
64namespace {
66 constexpr size_t iW = 0;
67 constexpr size_t iX = 1;
68 constexpr size_t iY = 2;
69 constexpr size_t iR2 = 3;
70
72 UncertainPerigeeCircle fitKarimaki(const double /*sw*/,
73 const Eigen::Matrix< double, 4, 1 >& a,
74 const Eigen::Matrix< double, 4, 4 >& c,
75 bool lineConstrained = false)
76 {
77 double q1, q2 = 0.0;
78 if (lineConstrained) {
79 q1 = c(iX, iY);
80 q2 = c(iX, iX) - c(iY, iY);
81 } else {
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);
84 }
85
86 double phi = 0.5 * atan2(2. * q1, q2);
87
88 double sinphi = sin(phi);
89 double cosphi = cos(phi);
90
91 double curv, I = 0.0;
92 if (lineConstrained) {
93 curv = 0.0; //line
94 I = sinphi * a(iX) - cosphi * a(iY);
95
96 } else {
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));
101
102 }
103
104 // Karimaki uses the opposite sign for phi in contrast to the convention of this framework !!!
105 phi += phi > 0 ? -M_PI : M_PI;
106 return UncertainPerigeeCircle(curv, phi, I);
107
108 }
109
110
111
112
113
115 double calcChi2Karimaki(const PerigeeCircle& parameters,
116 const double sw,
117 const Eigen::Matrix< double, 4, 4 >& c,
118 bool lineConstrained = false)
119 {
120 // Karimaki uses the opposite sign for phi in contrast to the convention of this framework !!!
121 const Vector2D vecPhi = -parameters.tangential();
122
123 const double cosphi = vecPhi.x();
124 const double sinphi = vecPhi.y();
125
126 if (lineConstrained) {
127 double chi2 = sw * (sinphi * sinphi * c(iX, iX) - 2. * sinphi * cosphi * c(iX, iY) + cosphi * cosphi * c(iY, iY));
128 return chi2;
129 } else {
130 // Terminology Karimaki used in the paper
131 const double rho = parameters.curvature();
132 const double d = parameters.impact();
133
134 const double u = 1 + d * rho;
135 const double kappa = 0.5 * rho / u;
136
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));
139 return chi2;
140 }
141
142 }
143
144
145
146 PerigeePrecision calcPrecisionKarimaki(const PerigeeCircle& parameters,
147 const Eigen::Matrix< double, 4, 4 >& s,
148 bool lineConstrained = false)
149 {
150 PerigeePrecision perigeePrecision;
151
152 const double curv = parameters.curvature();
153 const double I = parameters.impact();
154
155 // Karimaki uses the opposite sign for phi in contrast to the convention of this framework !!!
156 const Vector2D vecPhi = -parameters.tangential();
157
158 // Ternminology Karimaki using in the paper
159 const double cosphi = vecPhi.x();
160 const double sinphi = vecPhi.y();
161
162 const double ssphi = sinphi * sinphi;
163 const double scphi = sinphi * cosphi;
164 const double ccphi = cosphi * cosphi;
165
166 const double rho = curv;
167 const double d = I;
168
169 const double u = 1. + rho * d;
170
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.;
178
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);
183
184 } else {
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));
188
189 double sd = sinphi * s(iX, iR2) - cosphi * s(iY, iR2);
190
191 double saa = ssphi * s(iX, iX) - 2. * scphi * s(iX, iY) + ccphi * s(iY, iY);
192
193 // Not in the Karimaki paper, but factors a similar term.
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);
196
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;
201
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);
207 }
208 return perigeePrecision;
209 }
210}
211
212
213
215{
216 // Matrix of weighted sums
217 Eigen::Matrix< double, 4, 4> sNoL = getWXYRSumMatrix(observations2D);
218
219 // Matrix of averages
220 Eigen::Matrix<double, 4, 4> aNoL = sNoL / sNoL(iW);
221
222 // Measurement means
223 Eigen::Matrix<double, 4, 1> meansNoL = aNoL.row(iW);
224
225 // Covariance matrix
226 Eigen::Matrix<double, 4, 4> cNoL = aNoL - meansNoL * meansNoL.transpose();
227
228 // Determine NDF : Circle fit eats up to 3 degrees of freedom debpending on the constraints
229 size_t ndf = observations2D.size() - 2;
230
231 if (not isLineConstrained()) {
232 --ndf;
233 }
234
235 // Parameters to be fitted
236 UncertainPerigeeCircle resultCircle = fitKarimaki(sNoL(iW), meansNoL, cNoL, isLineConstrained());
237 double chi2 = calcChi2Karimaki(resultCircle, sNoL(iW), cNoL);
238
239 PerigeePrecision perigeePrecision = calcPrecisionKarimaki(resultCircle, sNoL, isLineConstrained());
240
241 // Use in pivotingin caset the matrix is not full rank as is for the constrained cases-
242 PerigeeCovariance perigeeCovariance = PerigeeUtil::covarianceFromPrecision(perigeePrecision);
243
244 resultCircle.setChi2(chi2);
245 resultCircle.setNDF(ndf);
246 resultCircle.setPerigeeCovariance(perigeeCovariance);
247 return resultCircle;
248}
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.
Definition: PerigeeCircle.h:36
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.
Definition: PlainMatrix.h:40
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...
Definition: Vector2D.h:32
double x() const
Getter for the x coordinate.
Definition: Vector2D.h:595
double y() const
Getter for the y coordinate.
Definition: Vector2D.h:605
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
static CovarianceMatrix covarianceFromPrecision(const PrecisionMatrix &prec)
Convert the precision matrix to the corresponding covariance matrix.