Belle II Software  release-08-01-10
RiemannsMethod Class Reference

Class implementing the Riemann fit for two dimensional trajectory circle. More...

#include <RiemannsMethod.h>

Public Member Functions

 RiemannsMethod ()
 Default constructor.
 
void update (CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
 Executes the fit and updates the trajectory parameters. This may render the information in the observation object.
 
bool isLineConstrained () const
 Getter for the indictor that lines should be fitted by this fitter.
 
bool isOriginConstrained () const
 Getter for the indictor that curves through the origin should be fitted by this fitter.
 
void setLineConstrained (bool constrained=true)
 Indicator if this fitter is setup to fit lines.
 
void setOriginConstrained (bool constrained=true)
 Indicator if this fitter is setup to fit curves through the origin.
 

Private Member Functions

void updateWithoutDriftLength (CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
 Executes the fit without using the drift length information. More...
 
void updateWithDriftLength (CDCTrajectory2D &trajectory2D, CDCObservations2D &observations2D) const
 Executes the fit using the drift length information. More...
 

Private Attributes

bool m_lineConstrained
 Memory for the flag indicating that lines should be fitted.
 
bool m_originConstrained
 Memory for the flag indicating that curves through the origin shall be fitter.
 

Detailed Description

Class implementing the Riemann fit for two dimensional trajectory circle.

Definition at line 20 of file RiemannsMethod.h.

Member Function Documentation

◆ updateWithDriftLength()

void updateWithDriftLength ( CDCTrajectory2D trajectory2D,
CDCObservations2D observations2D 
) const
private

Executes the fit using the drift length information.

This method is used if there is no drift length information is available from the observations.

Definition at line 144 of file RiemannsMethod.cc.

145 {
146  EigenObservationMatrix eigenObservation = getEigenObservationMatrix(&observations2D);
147  assert(eigenObservation.cols() == 4);
148  size_t nObservations = observations2D.size();
149 
150  //cout << "updateWithRightLeft : " << endl;
151  //observations always have the structure
152  /*
153  observables[0][0] == x of first point
154  observables[0][1] == y of first point
155  observables[0][2] == desired distance of first point
156  */
157 
158  Eigen::Matrix< double, Eigen::Dynamic, 1 > distances = eigenObservation.col(2);
159 
160  //cout << "distances : " << endl << distances << endl;
161 
162  if ((isLineConstrained()) && (isOriginConstrained())) {
163 
164 
165  Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 2);
166 
167  //all coordiates
168  //projectedPoints.col(0) = Eigen::Matrix<double,Eigen::Dynamic,1>::Constant(nObservations,1.0);
169  projectedPoints.col(0) = eigenObservation.col(0);
170  projectedPoints.col(1) = eigenObservation.col(1);
171  //projectedPoints.col(2) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
172 
173  Eigen::Matrix<double, 2, 1> parameters =
174  projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
175 
176  trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(0.0, parameters(0), parameters(1), 0.0)));
177  }
178 
179  else if ((! isLineConstrained()) && (isOriginConstrained())) {
180 
181  Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 3);
182 
183  //all coordiates
184  //projectedPoints.col(0) = Eigen::Matrix<double,Eigen::Dynamic,1>::Constant(1.0);
185  projectedPoints.col(0) = eigenObservation.col(0);
186  projectedPoints.col(1) = eigenObservation.col(1);
187  projectedPoints.col(2) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
188 
189  Eigen::Matrix<double, 3, 1> parameters =
190  projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
191 
192  trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(0.0, parameters(0), parameters(1), parameters(2))));
193  }
194 
195  else if ((isLineConstrained()) && (! isOriginConstrained())) {
196 
197  Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 3);
198 
199  //all coordiates
200  projectedPoints.col(0) = Eigen::Matrix<double, Eigen::Dynamic, 1>::Constant(nObservations, 1.0);
201  projectedPoints.col(1) = eigenObservation.col(0);
202  projectedPoints.col(2) = eigenObservation.col(1);
203  //projectedPoints.col(3) = eigenObservation.leftCols<2>().rowwise().squaredNorm()- distances.rowwise().squaredNorm();
204 
205  Eigen::Matrix<double, 3, 1> parameters =
206  projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
207 
208  trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(parameters(0), parameters(1), parameters(2), 0.0)));
209  //fit.setParameters(parameters(0),parameters(1),parameters(2),0.0);
210 
211  }
212 
213  else if ((! isLineConstrained()) && (! isOriginConstrained())) {
214 
215  Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > projectedPoints(nObservations, 4);
216 
217  //all coordiates
218  projectedPoints.col(0) = Eigen::Matrix<double, Eigen::Dynamic, 1>::Constant(nObservations, 1.0);
219  projectedPoints.col(1) = eigenObservation.col(0);
220  projectedPoints.col(2) = eigenObservation.col(1);
221  projectedPoints.col(3) = eigenObservation.leftCols<2>().rowwise().squaredNorm() - distances.rowwise().squaredNorm();
222 
223  //cout << "points : " << endl << projectedPoints << endl;
224 
225  Eigen::Matrix<double, 4, 1> parameters =
226  projectedPoints.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(distances);
227 
228  trajectory2D.setGlobalCircle(UncertainPerigeeCircle(PerigeeCircle::fromN(parameters(0), parameters(1), parameters(2),
229  parameters(3))));
230 
231  }
232 }
std::size_t size() const
Returns the number of observations stored.
void setGlobalCircle(const UncertainPerigeeCircle &perigeeCircle)
Setter for the generalized circle that describes the trajectory.
static PerigeeCircle fromN(double n0, double n1, double n2, double n3=0)
Constructor with the four parameters of the generalized circle.
bool isLineConstrained() const
Getter for the indictor that lines should be fitted by this fitter.
bool isOriginConstrained() const
Getter for the indictor that curves through the origin should be fitted by this fitter.
Adds an uncertainty matrix to the circle in perigee parameterisation.

◆ updateWithoutDriftLength()

void updateWithoutDriftLength ( CDCTrajectory2D trajectory2D,
CDCObservations2D observations2D 
) const
private

Executes the fit without using the drift length information.

This method is used if there is drift length information is available from the observations.

Definition at line 47 of file RiemannsMethod.cc.


The documentation for this class was generated from the following files: