1 #ifndef RecoTracker_TkSeedGenerator_FastCircleFit_h 2 #define RecoTracker_TkSeedGenerator_FastCircleFit_h 12 #define PRINT std::cout 15 #define PRINT LogTrace("FastCircleFit") 37 template <
typename P,
typename E>
39 const auto N = points.size();
52 std::array<float, N>
x;
53 std::array<float, N>
y;
54 std::array<float, N>
z;
55 std::array<float, N>
weight;
61 float x0()
const {
return x0_; }
62 float y0()
const {
return y0_; }
69 template <
typename P,
typename E,
typename C>
83 template <
typename P,
typename E,
typename C>
85 const auto N = points.size();
88 for (
size_t i = 0;
i <
N; ++
i) {
89 const auto&
point = points[
i];
90 const auto&
p =
point.basicVector();
99 const float w = 1.f / (
point.perp2() * phiErr2);
102 PRINT <<
" point " <<
p.x() <<
" " <<
p.y() <<
" transformed " <<
x[
i] <<
" " <<
y[
i] <<
" " <<
z[
i] <<
" weight " 105 const float invTotWeight = 1.f / std::accumulate(
weight.begin(),
weight.end(), 0.f);
106 PRINT <<
" invTotWeight " << invTotWeight;
109 for (
size_t i = 0;
i <
N; ++
i) {
115 mean *= invTotWeight;
119 for (
size_t i = 0;
i <
N; ++
i) {
120 const auto diff_x =
x[
i] -
mean[0];
121 const auto diff_y =
y[
i] -
mean[1];
122 const auto diff_z =
z[
i] -
mean[2];
127 A(0, 0) +=
w * diff_x * diff_x;
128 A(1, 0) +=
w * diff_y * diff_x;
129 A(1, 1) +=
w * diff_y * diff_y;
130 A(2, 0) +=
w * diff_z * diff_x;
131 A(2, 1) +=
w * diff_z * diff_y;
132 A(2, 2) +=
w * diff_z * diff_z;
133 PRINT <<
" i " <<
i <<
" A " <<
A(0, 0) <<
" " <<
A(0, 1) <<
" " <<
A(0, 2) << std::endl
134 <<
" " <<
A(1, 0) <<
" " <<
A(1, 1) <<
" " <<
A(1, 2) << std::endl
135 <<
" " <<
A(2, 0) <<
" " <<
A(2, 1) <<
" " <<
A(2, 2) << std::endl;
137 A *= 1. /
static_cast<float>(
N);
139 PRINT <<
" A " <<
A(0, 0) <<
" " <<
A(0, 1) <<
" " <<
A(0, 2) << std::endl
140 <<
" " <<
A(1, 0) <<
" " <<
A(1, 1) <<
" " <<
A(1, 2) << std::endl
141 <<
" " <<
A(2, 0) <<
" " <<
A(2, 1) <<
" " <<
A(2, 2) << std::endl;
144 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigen(
A);
145 const auto& eigenValues = eigen.eigenvalues();
146 const auto& eigenVectors = eigen.eigenvectors();
149 PRINT <<
" eigenvalues " << eigenValues[0] <<
" " << eigenValues[1] <<
" " << eigenValues[2] << std::endl;
151 PRINT <<
" eigen " << eigenVectors(0, 0) <<
" " << eigenVectors(0, 1) <<
" " << eigenVectors(0, 2) << std::endl
152 <<
" " << eigenVectors(1, 0) <<
" " << eigenVectors(1, 1) <<
" " << eigenVectors(1, 2) << std::endl
153 <<
" " << eigenVectors(2, 0) <<
" " << eigenVectors(2, 1) <<
" " << eigenVectors(2, 2) << std::endl;
156 auto n = eigenVectors.col(0);
157 PRINT <<
" n (unit) " <<
n[0] <<
" " <<
n[1] <<
" " <<
n[2] << std::endl;
161 PRINT <<
" c " <<
c << std::endl;
164 const auto tmp = 0.5f * 1.f /
n[2];
167 const float rho2 = (1 -
sqr(
n[2]) - 4 *
c *
n[2]) *
sqr(
tmp);
172 for (
size_t i = 0;
i <
N; ++
i) {
174 #if defined(MK_DEBUG) || defined(EDM_ML_DEBUG) 176 PRINT <<
" i " <<
i <<
" chi2 incr " << incr <<
" d(hit, circle) " <<
distance <<
" sigma "
FastCircleFit(const P &points, const E &errors)
FastCircleFit(const std::array< GlobalPoint, N > &points, const std::array< GlobalError, N > &errors)
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
std::pair< OmniClusterRef, TrackingParticleRef > P
#define declareDynArray(T, n, x)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point