Belle II Software development
RiemannsMethod.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/RiemannsMethod.h>
9
10#include <tracking/trackFindingCDC/fitting/EigenObservationMatrix.h>
11#include <tracking/trackFindingCDC/fitting/CDCObservations2D.h>
12
13#include <tracking/trackingUtilities/eventdata/trajectories/CDCTrajectory2D.h>
14
15#include <framework/logging/Logger.h>
16
17#include <Eigen/Eigen>
18#include <Eigen/Core>
19
20#include <cassert>
21
22using namespace Belle2;
23using namespace TrackFindingCDC;
24using namespace TrackingUtilities;
25
31
33 CDCObservations2D& observations2D) const
34{
35
36 if (observations2D.getNObservationsWithDriftRadius() > 0) {
37 updateWithDriftLength(trajectory2D, observations2D);
38 } else {
39 updateWithoutDriftLength(trajectory2D, observations2D);
40 }
41
42 EForwardBackward isCoaligned = observations2D.isCoaligned(trajectory2D);
43 if (isCoaligned == EForwardBackward::c_Backward) trajectory2D.reverse();
44}
45
46
47
49 CDCObservations2D& observations2D) const
50{
51 EigenObservationMatrix eigenObservation = getEigenObservationMatrix(&observations2D);
52 assert(eigenObservation.cols() == 4);
53 size_t nObservations = observations2D.size();
54
55 if (isLineConstrained()) {
56
57 //do a normal line fit
58 Eigen::Matrix<double, Eigen::Dynamic, 2> points = eigenObservation.leftCols<2>();
59
60 Eigen::Matrix<double, 1, 2> pointMean;
61 //RowVector2d pointMean;
62 // cppcheck-suppress constStatement
63 pointMean << 0.0, 0.0;
64 if (!(isOriginConstrained())) {
65 // subtract the offset from the origin
66 pointMean = points.colwise().mean();
67
68 points = points.rowwise() - pointMean;
69 }
70 Eigen::Matrix<double, 2, 2> covarianceMatrix = points.transpose() * points;
71
72 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 2, 2> > eigensolver(covarianceMatrix);
73
74 if (eigensolver.info() != Eigen::Success)
75 B2WARNING("SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
76
77 //the eigenvalues are generated in increasing order
78 //we are interested in the lowest one since we want to compute the normal vector of the plane
79
80 Eigen::Matrix<double, 2, 1> normalToLine = eigensolver.eigenvectors().col(0);
81
82 double offset = -pointMean * normalToLine;
83
84 // set the generalized circle parameters
85 // last set to zero constrains to a line
86 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(offset, normalToLine(0), normalToLine(1), 0)));
87
88 } else {
89
90 //lift the points to the projection space
91 Eigen::Matrix<double, Eigen::Dynamic, 3> projectedPoints(nObservations, 3);
92
93 projectedPoints.col(0) = eigenObservation.col(0);
94 projectedPoints.col(1) = eigenObservation.col(1);
95 projectedPoints.col(2) = eigenObservation.leftCols<2>().rowwise().squaredNorm();
96
97
98 Eigen::Matrix<double, 1, 3> pointMean;
99 // cppcheck-suppress constStatement
100 pointMean << 0.0, 0.0, 0.0;
101 if (!(isOriginConstrained())) {
102 // subtract the offset from the origin
103 pointMean = projectedPoints.colwise().mean();
104
105 projectedPoints = projectedPoints.rowwise() - pointMean;
106 }
107
108 Eigen::Matrix<double, 3, 3> covarianceMatrix = projectedPoints.transpose() * projectedPoints;
109
110 Eigen::SelfAdjointEigenSolver< Eigen::Matrix<double, 3, 3> > eigensolver(covarianceMatrix);
111
112 if (eigensolver.info() != Eigen::Success)
113 B2WARNING("Eigen::SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
114
115 //the eigenvalues are generated in increasing order
116 //we are interested in the lowest one since we want to compute the normal vector of the plane
117
118 Eigen::Matrix<double, 3, 1> normalToPlane = eigensolver.eigenvectors().col(0);
119
120 double offset = -pointMean * normalToPlane;
121
122 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(offset, normalToPlane(0), normalToPlane(1),
123 normalToPlane(2))));
124 //fit.setParameters();
125
126 }
127
128 //check if the orientation is alright
129 Vector2D directionAtCenter = trajectory2D.getFlightDirection2D(Vector2D(0.0, 0.0));
130
131
132 size_t voteForChangeSign = 0;
133 for (size_t iPoint = 0; iPoint < nObservations; ++iPoint) {
134 double pointInSameDirection = eigenObservation(iPoint, 0) * directionAtCenter.x() +
135 eigenObservation(iPoint, 1) * directionAtCenter.y();
136 if (pointInSameDirection < 0) ++voteForChangeSign;
137 }
138
139 if (voteForChangeSign > nObservations / 2.0) trajectory2D.reverse();
140
141}
142
143
144
146{
147 EigenObservationMatrix eigenObservation = getEigenObservationMatrix(&observations2D);
148 assert(eigenObservation.cols() == 4);
149 size_t nObservations = observations2D.size();
150
151 //cout << "updateWithRightLeft : " << endl;
152 //observations always have the structure
153 /*
154 observables[0][0] == x of first point
155 observables[0][1] == y of first point
156 observables[0][2] == desired distance of first point
157 */
158
159 Eigen::Matrix< double, Eigen::Dynamic, 1 > distances = eigenObservation.col(2);
160
161 //cout << "distances : " << endl << distances << endl;
162
163 if ((isLineConstrained()) && (isOriginConstrained())) {
164
165
166 Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 2);
167
168 //all coordinates
169 //projectedPoints.col(0) = Eigen::Matrix<double,Eigen::Dynamic,1>::Constant(nObservations,1.0);
170 projectedPoints.col(0) = eigenObservation.col(0);
171 projectedPoints.col(1) = eigenObservation.col(1);
172 //projectedPoints.col(2) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
173
174 Eigen::Matrix<double, 2, 1> parameters =
175 projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
176
177 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(0.0, parameters(0), parameters(1), 0.0)));
178 }
179
180 else if ((! isLineConstrained()) && (isOriginConstrained())) {
181
182 Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 3);
183
184 //all coordinates
185 //projectedPoints.col(0) = Eigen::Matrix<double,Eigen::Dynamic,1>::Constant(1.0);
186 projectedPoints.col(0) = eigenObservation.col(0);
187 projectedPoints.col(1) = eigenObservation.col(1);
188 projectedPoints.col(2) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
189
190 Eigen::Matrix<double, 3, 1> parameters =
191 projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
192
193 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(0.0, parameters(0), parameters(1), parameters(2))));
194 }
195
196 else if ((isLineConstrained()) && (! isOriginConstrained())) {
197
198 Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 3);
199
200 //all coordinates
201 projectedPoints.col(0) = Eigen::Matrix<double, Eigen::Dynamic, 1>::Constant(nObservations, 1.0);
202 projectedPoints.col(1) = eigenObservation.col(0);
203 projectedPoints.col(2) = eigenObservation.col(1);
204 //projectedPoints.col(3) = eigenObservation.leftCols<2>().rowwise().squaredNorm()- distances.rowwise().squaredNorm();
205
206 Eigen::Matrix<double, 3, 1> parameters =
207 projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
208
209 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(parameters(0), parameters(1), parameters(2), 0.0)));
210 //fit.setParameters(parameters(0),parameters(1),parameters(2),0.0);
211
212 }
213
214 else if ((! isLineConstrained()) && (! isOriginConstrained())) {
215
216 Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 4);
217
218 //all coordinates
219 projectedPoints.col(0) = Eigen::Matrix<double, Eigen::Dynamic, 1>::Constant(nObservations, 1.0);
220 projectedPoints.col(1) = eigenObservation.col(0);
221 projectedPoints.col(2) = eigenObservation.col(1);
222 projectedPoints.col(3) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
223
224 //cout << "points : " << endl << projectedPoints << endl;
225
226 Eigen::Matrix<double, 4, 1> parameters =
227 projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
228
229 trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(parameters(0), parameters(1), parameters(2),
230 parameters(3))));
231
232 }
233}
Class serving as a storage of observed drift circles to present to the Riemann fitter.
std::size_t size() const
Returns the number of observations stored.
TrackingUtilities::EForwardBackward isCoaligned(const TrackingUtilities::CDCTrajectory2D &trajectory2D) const
Checks if the last observation in the vector lies at greater or lower travel distance than the last o...
std::size_t getNObservationsWithDriftRadius() const
Returns the number of observations having a drift radius radius.
void update(TrackingUtilities::CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
Executes the fit and updates the trajectory parameters. This may render the information in the observ...
void updateWithoutDriftLength(TrackingUtilities::CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
Executes the fit without using the drift length information.
bool isLineConstrained() const
Getter for the indicator that lines should be fitted by this fitter.
bool m_originConstrained
Memory for the flag indicating that curves through the origin shall be fitter.
void updateWithDriftLength(TrackingUtilities::CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
Executes the fit using the drift length information.
bool m_lineConstrained
Memory for the flag indicating that lines should be fitted.
bool isOriginConstrained() const
Getter for the indicator that curves through the origin should be fitted by this fitter.
Particle trajectory as it is seen in xy projection represented as a circle.
void reverse()
Reverses the trajectory in place.
void setGlobalCircle(const UncertainPerigeeCircle &perigeeCircle)
Setter for the generalized circle that describes the trajectory.
Vector2D getFlightDirection2D(const Vector2D &point) const
Get the unit direction of flight at the given point, where arcLength2D = 0.
static PerigeeCircle fromN(double n0, double n1, double n2, double n3=0)
Constructor with the four parameters of the generalized circle.
Adds an uncertainty matrix to the circle in perigee parameterisation.
A two dimensional vector which is equipped with functions for correct handling of orientation relate...
Definition Vector2D.h:36
double x() const
Getter for the x coordinate.
Definition Vector2D.h:626
double y() const
Getter for the y coordinate.
Definition Vector2D.h:641
Abstract base class for different kinds of events.