Belle II Software development
CDCSZFitter.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/CDCSZFitter.h>
9
10#include <tracking/trackFindingCDC/fitting/EigenObservationMatrix.h>
11#include <tracking/trackFindingCDC/fitting/CDCSZObservations.h>
12#include <tracking/trackFindingCDC/fitting/CDCObservations2D.h>
13
14#include <tracking/trackFindingCDC/eventdata/tracks/CDCTrack.h>
15#include <tracking/trackFindingCDC/eventdata/tracks/CDCSegmentPair.h>
16
17#include <tracking/trackFindingCDC/eventdata/segments/CDCSegment3D.h>
18#include <tracking/trackFindingCDC/eventdata/segments/CDCSegment2D.h>
19
20#include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectorySZ.h>
21#include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectory2D.h>
22
23#include <tracking/trackFindingCDC/geometry/SZParameters.h>
24
25#include <tracking/trackFindingCDC/numerics/EigenView.h>
26
27#include <framework/logging/Logger.h>
28
29#include <Eigen/Eigen>
30#include <Eigen/QR>
31#include <Eigen/Core>
32
33using namespace Belle2;
34using namespace TrackFindingCDC;
35
37{
38 static CDCSZFitter szFitter;
39 return szFitter;
40}
41
42namespace {
43 UncertainSZLine fitZ(const Eigen::Matrix<double, 3, 3>& sumMatrixWSZ)
44 {
45 // Solve the normal equation X * n = y
46 Eigen::Matrix<double, 2, 2> sumMatrixWS = sumMatrixWSZ.block<2, 2>(0, 0);
47 Eigen::Matrix<double, 2, 1> sumVectorZOverWS = sumMatrixWSZ.block<2, 1>(0, 2);
48 Eigen::Matrix<double, 2, 1> nZOverWS = sumMatrixWS.llt().solve(sumVectorZOverWS);
49 double chi2 = sumMatrixWSZ(2, 2) - nZOverWS.transpose() * sumVectorZOverWS;
50
51 using namespace NSZParameterIndices;
52 SZParameters szParameters;
53 szParameters(c_TanL) = nZOverWS(1);
54 szParameters(c_Z0) = nZOverWS(0);
55
56 SZPrecision szPrecision;
57 szPrecision(c_TanL, c_TanL) = sumMatrixWS(1, 1);
58 szPrecision(c_Z0, c_TanL) = sumMatrixWS(0, 1);
59 szPrecision(c_TanL, c_Z0) = sumMatrixWS(1, 0);
60 szPrecision(c_Z0, c_Z0) = sumMatrixWS(0, 0);
61
62 SZCovariance szCovariance = SZUtil::covarianceFromFullPrecision(szPrecision);
63 return UncertainSZLine(szParameters, szCovariance, chi2);
64 }
65
66 // Declare function as currently unused to avoid compiler warning
67 UncertainSZLine fitSZ(const Eigen::Matrix<double, 3, 3>& sumMatrixWSZ) __attribute__((__unused__));
68
70 UncertainSZLine fitSZ(const Eigen::Matrix<double, 3, 3>& sumMatrixWSZ)
71 {
72
73 // Matrix of averages
74 Eigen::Matrix<double, 3, 3> averageMatrixWSZ = sumMatrixWSZ / sumMatrixWSZ(0);
75
76 // Measurement means
77 Eigen::Matrix<double, 3, 1> meansWSZ = averageMatrixWSZ.row(0);
78
79 // Covariance matrix
80 Eigen::Matrix<double, 3, 3> covMatrixWSZ = averageMatrixWSZ - meansWSZ * meansWSZ.transpose();
81
82 Eigen::Matrix<double, 2, 2> covMatrixSZ = covMatrixWSZ.block<2, 2>(1, 1);
83 Eigen::Matrix<double, 2, 1> meansSZ = meansWSZ.segment<2>(1);
84
85 Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 2, 2>> eigensolver(covMatrixSZ);
86 if (eigensolver.info() != Eigen::Success) {
87 B2WARNING(
88 "SelfAdjointEigenSolver could not compute the eigen values of the observation matrix");
89 }
90
91 // the eigenvalues are generated in increasing order
92 // we are interested in the lowest one since we want to compute the normal vector of the line
93 // here
94 Eigen::Matrix<double, 2, 1> nSZ = eigensolver.eigenvectors().col(0);
95 B2INFO("nSZ " << nSZ);
96 nSZ /= std::copysign(nSZ.norm(), -nSZ(1)); // Making n2 negative to normalize to forward along s
97 Eigen::Matrix<double, 3, 1> nWSZ;
98 // cppcheck-suppress constStatement
99 nWSZ << -meansSZ.transpose() * nSZ, nSZ;
100 B2INFO("nWSZ " << nWSZ);
101
102 double chi2 = nWSZ.transpose() * sumMatrixWSZ * nWSZ;
103
104 Eigen::Matrix<double, 3, 3> precN = sumMatrixWSZ;
105
106 using namespace NSZParameterIndices;
107 SZParameters szParameters;
108 szParameters(c_TanL) = -nWSZ(1) / nWSZ(2);
109 szParameters(c_Z0) = -nWSZ(0) / nWSZ(2);
110
111 Eigen::Matrix<double, 3, 2> ambiguitySZToN = Eigen::Matrix<double, 3, 2>::Zero();
112 ambiguitySZToN(0, c_Z0) = -nWSZ(2);
113 ambiguitySZToN(0, c_TanL) = -nWSZ(1) * nWSZ(0);
114 ambiguitySZToN(1, c_TanL) = -nWSZ(2) * nWSZ(2) * nWSZ(2);
115 ambiguitySZToN(2, c_TanL) = nWSZ(1) * nWSZ(2) * nWSZ(2);
116
117 SZPrecision szPrecision;
118 mapToEigen(szPrecision) = ambiguitySZToN.transpose() * precN * ambiguitySZToN;
119
120 SZCovariance szCovariance = SZUtil::covarianceFromFullPrecision(szPrecision);
121 return UncertainSZLine(szParameters, szCovariance, chi2);
122 }
123}
124
126{
127 const bool onlyStereo = true;
128 CDCSZObservations observationsSZ(EFitVariance::c_Proper, onlyStereo);
129 observationsSZ.appendRange(track);
130 if (observationsSZ.size() > 3) {
131 CDCTrajectorySZ szTrajectory;
132 update(szTrajectory, observationsSZ);
133 return szTrajectory;
134 } else {
136 }
137}
138
140 const CDCTrajectory2D& axialTrajectory2D) const
141{
142 B2ASSERT("Expected stereo segment", not stereoSegment.isAxial());
143
144 CDCTrajectorySZ trajectorySZ;
145 update(trajectorySZ, stereoSegment, axialTrajectory2D);
146 return trajectorySZ;
147}
148
150{
151 CDCSZObservations observationsSZ;
152 observationsSZ.appendRange(segment3D);
153 return fit(std::move(observationsSZ));
154}
155
157{
158 CDCTrajectorySZ trajectorySZ;
159 update(trajectorySZ, observationsSZ);
160 return trajectorySZ;
161}
162
163void CDCSZFitter::update(const CDCSegmentPair& segmentPair) const
164{
165 const CDCSegment2D* ptrStereoSegment = segmentPair.getStereoSegment();
166 const CDCSegment2D* ptrAxialSegment = segmentPair.getAxialSegment();
167
168 assert(ptrStereoSegment);
169 assert(ptrAxialSegment);
170
171 const CDCSegment2D& stereoSegment = *ptrStereoSegment;
172 const CDCSegment2D& axialSegment = *ptrAxialSegment;
173 const CDCTrajectory2D& axialTrajectory2D = axialSegment.getTrajectory2D();
174
175 CDCTrajectorySZ trajectorySZ;
176 update(trajectorySZ, stereoSegment, axialTrajectory2D);
177
178 CDCTrajectory3D trajectory3D(axialTrajectory2D, trajectorySZ);
179 segmentPair.setTrajectory3D(trajectory3D);
180}
181
183{
184 CDCSZObservations szObservations;
185 for (size_t i = 0; i < observations2D.size(); ++i) {
186 const double s = observations2D.getX(i);
187 const double z = observations2D.getY(i);
188 szObservations.fill(s, z);
189 }
190 return fit(std::move(szObservations));
191}
192
194 const CDCSegment2D& stereoSegment,
195 const CDCTrajectory2D& axialTrajectory2D) const
196{
197 B2ASSERT("Expected stereo segment", not stereoSegment.isAxial());
198
199 // recostruct the stereo segment
200 CDCSZObservations observationsSZ;
201 for (const CDCRecoHit2D& recoHit2D : stereoSegment) {
202 CDCRecoHit3D recoHit3D = CDCRecoHit3D::reconstruct(recoHit2D, axialTrajectory2D);
203 observationsSZ.append(recoHit3D);
204 }
205
206 update(trajectorySZ, observationsSZ);
207}
208
209void CDCSZFitter::update(CDCTrajectorySZ& trajectorySZ, CDCSZObservations& observationsSZ) const
210{
211 trajectorySZ.clear();
212 if (observationsSZ.size() < 3) {
213 return;
214 }
215
216 // Determine NDF : Line fit eats up 2 degrees of freedom.
217 size_t ndf = observationsSZ.size() - 2;
218
219 // Matrix of weighted sums
220 Eigen::Matrix<double, 3, 3> sumMatrixWSZ = getWSZSumMatrix(observationsSZ);
221 UncertainSZLine uncertainSZLine = fitZ(sumMatrixWSZ);
222
223 uncertainSZLine.setNDF(ndf);
224 trajectorySZ.setSZLine(uncertainSZLine);
225}
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.
std::size_t size() const
Returns the number of observations stored.
Class representing a two dimensional reconstructed hit in the central drift chamber.
Definition: CDCRecoHit2D.h:47
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:52
static CDCRecoHit3D reconstruct(const CDCRecoHit2D &recoHit2D, const CDCTrajectory2D &trajectory2D)
Reconstructs the three dimensional hit from the two dimensional and the two dimensional trajectory.
Definition: CDCRecoHit3D.cc:56
Class implementing the z coordinate over travel distance line fit.
Definition: CDCSZFitter.h:27
void update(const CDCSegmentPair &segmentPair) const
Updates the trajectory of the axial stereo segment pair inplace.
Definition: CDCSZFitter.cc:163
static const CDCSZFitter & getFitter()
Getter for a standard sz line fitter instance.
Definition: CDCSZFitter.cc:36
CDCTrajectorySZ fitWithStereoHits(const CDCTrack &track) const
Returns the fitted sz trajectory of the track with the z-information of all stereo hits of the number...
Definition: CDCSZFitter.cc:125
CDCTrajectorySZ fit(const CDCSegment2D &stereoSegment, const CDCTrajectory2D &axialTrajectory2D) const
Returns a fitted trajectory.
Definition: CDCSZFitter.cc:139
Class serving as a storage of observed sz positions to present to the sz line fitters.
std::size_t appendRange(const std::vector< CDCRecoHit3D > &recoHit3Ds)
Appends all reconstructed hits from the three dimensional track.
std::size_t fill(double s, double z, double weight=1.0)
Appends the observed position.
std::size_t append(const CDCRecoHit3D &recoHit3D)
Appends the observed position.
std::size_t size() const
Returns the number of observations stored.
A reconstructed sequence of two dimensional hits in one super layer.
Definition: CDCSegment2D.h:39
A segment consisting of three dimensional reconstructed hits.
Definition: CDCSegment3D.h:26
Class representing a pair of one reconstructed axial segement and one stereo segment in adjacent supe...
const CDCSegment2D * getAxialSegment() const
Getter for the axial segment.
void setTrajectory3D(const CDCTrajectory3D &trajectory3D) const
Setter for the three dimensional trajectory.
const CDCSegment2D * getStereoSegment() const
Getter for the stereo segment.
bool isAxial() const
Indicator if the underlying wires are axial.
Definition: CDCSegment.h:45
CDCTrajectory2D & getTrajectory2D() const
Getter for the two dimensional trajectory fitted to the segment.
Definition: CDCSegment.h:69
Class representing a sequence of three dimensional reconstructed hits.
Definition: CDCTrack.h:41
Particle trajectory as it is seen in xy projection represented as a circle.
Particle full three dimensional trajectory.
Linear trajectory in sz space.
static CDCTrajectorySZ basicAssumption()
Constucts a basic assumption, what the z0 start position and the sz slope are, including some broad v...
void setSZLine(const UncertainSZLine &szLine)
Setter for the line in sz space.
void clear()
Clears all information from this trajectory line.
A matrix implementation to be used as an interface typ through out the track finder.
Definition: PlainMatrix.h:40
A line in sz where s is the transverse travel distance as seen in the xy projection with uncertaintie...
void setNDF(std::size_t ndf)
Setter for the number of degrees of freediom used in the line fit.
Abstract base class for different kinds of events.
static CovarianceMatrix covarianceFromFullPrecision(const PrecisionMatrix &prec)
Convert the precision matrix to the corresponding covariance matrix.