Belle II Software development
ExtendedRiemannsMethod.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/ExtendedRiemannsMethod.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 <tracking/trackFindingCDC/numerics/EigenView.h>
20
21#include <framework/logging/Logger.h>
22
23#include <Eigen/Eigen>
24#include <Eigen/Core>
25
26using namespace Belle2::TrackFindingCDC;
27
29 : m_lineConstrained(false)
30 , m_originConstrained(false)
31{
32}
33
35 CDCObservations2D& observations2D) const
36{
37 size_t nObservations = observations2D.size();
38 trajectory2D.clear();
39 if (not nObservations) return;
40
41 Vector2D origin = Vector2D(0.0, 0.0);
42 Vector2D centralPoint = observations2D.getCentralPoint();
43
44 const Vector2D& ref = isOriginConstrained() ? origin : centralPoint;
45 observations2D.passiveMoveBy(ref);
46
47 UncertainPerigeeCircle perigeeCircle = fitInternal(observations2D);
48
49 double frontX = observations2D.getX(0);
50 double frontY = observations2D.getY(0);
51 Vector2D frontPos(frontX, frontY);
52
53 double backX = observations2D.getX(nObservations - 1);
54 double backY = observations2D.getY(nObservations - 1);
55 Vector2D backPos(backX, backY);
56
57 Vector2D overPos(0, 0);
58 double totalPerps = (perigeeCircle->arcLengthBetween(frontPos, overPos) +
59 perigeeCircle->arcLengthBetween(overPos, backPos));
60
61 if (totalPerps < 0) {
62 perigeeCircle.reverse();
63 }
64
65 trajectory2D = CDCTrajectory2D(ref, perigeeCircle);
66}
67
68namespace {
69
73 enum EParabolicIndices {
74 iW = 0,
75 iX = 1,
76 iY = 2,
77 iR2 = 3,
78 iL = 4
79 };
80 }
82 PerigeeCircle fit(const Eigen::Matrix< double, 5, 5 >& sumMatrix,
83 bool lineConstrained = false,
84 bool originConstrained = false)
85 {
86 using namespace NParabolicParameterIndices;
87 // Solve the normal equation X * n = y
88 if (lineConstrained) {
89 if (originConstrained) {
90 Eigen::Matrix< double, 2, 2> X = sumMatrix.block<2, 2>(1, 1);
91 Eigen::Matrix< double, 2, 1> y = sumMatrix.block<2, 1>(1, iL);
92 Eigen::Matrix< double, 2, 1> n = X.ldlt().solve(y);
93 return PerigeeCircle::fromN(0.0, n(0), n(1), 0.0);
94
95 } else {
96 Eigen::Matrix< double, 3, 3> X = sumMatrix.block<3, 3>(0, 0);
97 Eigen::Matrix< double, 3, 1> y = sumMatrix.block<3, 1>(0, iL);
98 Eigen::Matrix< double, 3, 1> n = X.ldlt().solve(y);
99 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), 0.0);
100
101 }
102
103 } else {
104 if (originConstrained) {
105 Eigen::Matrix< double, 3, 3> X = sumMatrix.block<3, 3>(1, 1);
106 Eigen::Matrix< double, 3, 1> y = sumMatrix.block<3, 1>(1, iL);
107 Eigen::Matrix< double, 3, 1> n = X.ldlt().solve(y);
108 return PerigeeCircle::fromN(0.0, n(0), n(1), n(2));
109
110 } else {
111 Eigen::Matrix< double, 4, 4> X = sumMatrix.block<4, 4>(0, 0);
112 Eigen::Matrix< double, 4, 1> y = sumMatrix.block<4, 1>(0, iL);
113 Eigen::Matrix< double, 4, 1> n = X.ldlt().solve(y);
114 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), n(iR2));
115
116 }
117 }
118 }
119
120
122 PerigeeCircle fit(Eigen::Matrix< double, 4, 4 > sumMatrix,
123 bool lineConstrained = false,
124 bool originConstrained = false)
125 {
126 // Solve the normal equation min_n n^T * X * n
127 // n is the smallest eigenvector
128
129 using namespace NParabolicParameterIndices;
130 if (lineConstrained) {
131 if (originConstrained) {
132 Eigen::Matrix< double, 2, 2> X = sumMatrix.block<2, 2>(1, 1);
133 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 2, 2> > eigensolver(X);
134 Eigen::Matrix<double, 2, 1> n = eigensolver.eigenvectors().col(0);
135 if (eigensolver.info() != Eigen::Success) {
136 B2WARNING("SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
137 }
138 return PerigeeCircle::fromN(0.0, n(0), n(1), 0.0);
139
140 } else {
141 Eigen::Matrix< double, 3, 3> X = sumMatrix.block<3, 3>(0, 0);
142 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 3, 3> > eigensolver(X);
143 Eigen::Matrix<double, 3, 1> n = eigensolver.eigenvectors().col(0);
144 if (eigensolver.info() != Eigen::Success) {
145 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
146 }
147 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), 0.0);
148
149 }
150 } else {
151 if (originConstrained) {
152 Eigen::Matrix< double, 3, 3> X = sumMatrix.block<3, 3>(1, 1);
153 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 3, 3> > eigensolver(X);
154 Eigen::Matrix<double, 3, 1> n = eigensolver.eigenvectors().col(0);
155 if (eigensolver.info() != Eigen::Success) {
156 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
157 }
158 return PerigeeCircle::fromN(0.0, n(0), n(1), n(2));
159
160 } else {
161 Eigen::Matrix< double, 4, 4 > X = sumMatrix.block<4, 4>(0, 0);
162 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 4, 4> > eigensolver(X);
163 Eigen::Matrix<double, 4, 1> n = eigensolver.eigenvectors().col(0);
164 if (eigensolver.info() != Eigen::Success) {
165 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
166 }
167
168 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), n(iR2));
169 }
170 }
171 }
172
174 PerigeeCircle fitSeperateOffset(Eigen::Matrix< double, 4, 1 > means,
175 Eigen::Matrix< double, 4, 4 > c,
176 bool lineConstrained = false)
177 {
178 // Solve the normal equation min_n n^T * c * n
179 // for the plane normal and move the plain by the offset
180 // n is the smallest eigenvector
181 using namespace NParabolicParameterIndices;
182 if (lineConstrained) {
183 Eigen::Matrix< double, 2, 2> X = c.block<2, 2>(1, 1);
184 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 2, 2> > eigensolver(X);
185 if (eigensolver.info() != Eigen::Success) {
186 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
187 }
188
189 //the eigenvalues are generated in increasing order
190 //we are interested in the lowest one since we want to compute the normal vector of the plane
191 Eigen::Matrix<double, 4, 1> n;
192 n.middleRows<2>(iX) = eigensolver.eigenvectors().col(0);
193 n(iW) = -means.middleRows<2>(iX).transpose() * n.middleRows<2>(iX);
194 n(iR2) = 0.;
195 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), n(iR2));
196
197 } else {
198 Eigen::Matrix< double, 3, 3> X = c.block<3, 3>(1, 1);
199 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 3, 3> > eigensolver(X);
200 if (eigensolver.info() != Eigen::Success) {
201 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
202 }
203
204 // the eigenvalues are generated in increasing order
205 // we are interested in the lowest one since we want to compute the normal vector of the plane
206 Eigen::Matrix<double, 4, 1> n;
207 n.middleRows<3>(iX) = eigensolver.eigenvectors().col(0);
208 n(iW) = -means.middleRows<3>(iX).transpose() * n.middleRows<3>(iX);
209 return PerigeeCircle::fromN(n(iW), n(iX), n(iY), n(iR2));
210
211 }
212 }
213
214
216 double calcChi2(const PerigeeCircle& parameters,
217 const Eigen::Matrix< double, 5, 5 >& s)
218 {
219 Eigen::Matrix<double, 5, 1> n;
220 n(0) = parameters.n0();
221 n(1) = parameters.n1();
222 n(2) = parameters.n2();
223 n(3) = parameters.n3();
224 n(4) = -1;
225
226 double chi2 = n.transpose() * s * n;
227 return chi2;
228 }
229
230
232 double calcChi2(const PerigeeCircle& parameters,
233 const Eigen::Matrix< double, 4, 4 >& s)
234 {
235 Eigen::Matrix<double, 4, 1> n;
236 n(0) = parameters.n0();
237 n(1) = parameters.n1();
238 n(2) = parameters.n2();
239 n(3) = parameters.n3();
240
241 double chi2 = n.transpose() * s * n;
242 return chi2;
243 }
244
245
246 PerigeePrecision calcPrecision(const PerigeeCircle& parameters,
247 const Eigen::Matrix< double, 4, 4 >& s,
248 bool lineConstrained = false,
249 bool originConstrained = false)
250 {
251 const double impact = parameters.impact();
252 const Vector2D& phi0Vec = parameters.tangential();
253 const double curvature = parameters.curvature();
254
255 using namespace NPerigeeParameterIndices;
256 Eigen::Matrix<double, 4, 3> ambiguity = Eigen::Matrix<double, 4, 3> ::Zero();
257
258 ambiguity(0, c_Curv) = impact * impact / 2;
259 ambiguity(1, c_Curv) = phi0Vec.y() * impact;
260 ambiguity(2, c_Curv) = -phi0Vec.x() * impact;
261 ambiguity(3, c_Curv) = 1.0 / 2.0;
262
263 ambiguity(0, c_Phi0) = 0;
264 ambiguity(1, c_Phi0) = phi0Vec.x() * (1 + curvature * impact);
265 ambiguity(2, c_Phi0) = phi0Vec.y() * (1 + curvature * impact);
266 ambiguity(3, c_Phi0) = 0;
267
268 ambiguity(0, c_I) = 1 + curvature * impact;
269 ambiguity(1, c_I) = phi0Vec.y() * curvature;
270 ambiguity(2, c_I) = -phi0Vec.x() * curvature;
271 ambiguity(3, c_I) = 0;
272
273 Eigen::Matrix<double, 3, 3> perigeePrecision = ambiguity.transpose() * s * ambiguity;
274 // Zero out the unfitted parameters from the precision matrix
275 if (lineConstrained) {
276 perigeePrecision.row(c_Curv) = Eigen::Matrix<double, 1, 3>::Zero();
277 perigeePrecision.col(c_Curv) = Eigen::Matrix<double, 3, 1>::Zero();
278 }
279
280 if (originConstrained) {
281 perigeePrecision.row(c_I) = Eigen::Matrix<double, 1, 3>::Zero();
282 perigeePrecision.col(c_I) = Eigen::Matrix<double, 3, 1>::Zero();
283 }
284
285 PerigeePrecision result;
286 mapToEigen(result) = perigeePrecision;
287 return result;
288 }
289
290}
291
292
293
295{
296 using namespace NParabolicParameterIndices;
297
298 // Matrix of weighted sums
299 Eigen::Matrix< double, 5, 5 > s = getWXYRLSumMatrix(observations2D);
300
301 // The same as above without drift lengths
302 Eigen::Matrix<double, 4, 4> sNoL = s.block<4, 4>(0, 0);
303
304 // Determine NDF : Circle fit eats up to 3 degrees of freedom depending on the constraints
305 size_t ndf = observations2D.size() - 1;
306
307 if (not isOriginConstrained()) {
308 --ndf;
309 }
310
311 if (not isLineConstrained()) {
312 --ndf;
313 }
314
315 // Parameters to be fitted
316 UncertainPerigeeCircle resultCircle;
317 double chi2 = 0;
318
319 size_t nObservationsWithDriftRadius = observations2D.getNObservationsWithDriftRadius();
320 if (nObservationsWithDriftRadius > 0) {
321 resultCircle = UncertainPerigeeCircle(::fit(s, isLineConstrained(), isOriginConstrained()));
322 chi2 = calcChi2(resultCircle, s);
323 } else {
324 if (not isOriginConstrained()) {
325 // Alternative implementation for comparison
326
327 // Matrix of averages
328 Eigen::Matrix< double, 4, 4> aNoL = sNoL / sNoL(iW);
329
330 // Measurement means
331 Eigen::Matrix< double, 4, 1> meansNoL = aNoL.row(iW);
332
333 // Covariance matrix
334 Eigen::Matrix< double, 4, 4> cNoL = aNoL - meansNoL * meansNoL.transpose();
335
336 resultCircle = UncertainPerigeeCircle(fitSeperateOffset(meansNoL, cNoL, isLineConstrained()));
337
338 } else {
339 resultCircle = UncertainPerigeeCircle(::fit(sNoL, isLineConstrained(), isOriginConstrained()));
340 }
341
342 chi2 = calcChi2(resultCircle, sNoL);
343 }
344
345 // Covariance calculation does not need the drift lengths, which is why we do not forward them.
346 PerigeePrecision perigeePrecision =
347 calcPrecision(resultCircle, sNoL, isLineConstrained(), isOriginConstrained());
348
349 // Use in pivoting in case the matrix is not full rank as it is for the constrained cases
350 PerigeeCovariance perigeeCovariance = PerigeeUtil::covarianceFromPrecision(perigeePrecision);
351
352 resultCircle.setNDF(ndf);
353 resultCircle.setChi2(chi2);
354 resultCircle.setPerigeeCovariance(perigeeCovariance);
355 return resultCircle;
356}
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.
std::size_t getNObservationsWithDriftRadius() const
Returns the number of observations having a drift radius radius.
Particle trajectory as it is seen in xy projection represented as a circle.
void clear()
Clears all information from this trajectory.
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.
ExtendedRiemannsMethod()
Constructor setting the default constraints.
void update(CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
Executes the fit and updates the trajectory parameters This may render the information in the observa...
bool isOriginConstrained() const
Getter for the indicator that curves through the origin should be fitted by this fitter.
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.
static PerigeeCircle fromN(double n0, double n1, double n2, double n3=0)
Constructor with the four parameters of the generalized 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
Namespace to hide the contained enum constants.
static CovarianceMatrix covarianceFromPrecision(const PrecisionMatrix &prec)
Convert the precision matrix to the corresponding covariance matrix.