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;
56 calculate(points, errors, x, y, z, weight);
61 float x0()
const {
return x0_; }
62 float y0()
const {
return y0_; }
69 template <
typename P,
typename E,
typename C>
81 template <
typename P,
typename E,
typename C>
83 const auto N = points.size();
86 for(
size_t i=0;
i<
N; ++
i) {
87 const auto&
point = points[
i];
88 const auto&
p =
point.basicVector();
96 const float phiErr2 = errors[
i].phierr(
point);
97 const float w = 1.f/(
point.perp2()*phiErr2);
100 PRINT <<
" point " <<
p.x()
102 <<
" transformed " << x[
i]
108 const float invTotWeight = 1.f/std::accumulate(weight.begin(), weight.end(), 0.f);
109 PRINT <<
" invTotWeight " << invTotWeight;
111 Eigen::Vector3f
mean = Eigen::Vector3f::Zero();
112 for(
size_t i=0;
i<
N; ++
i) {
113 const float w = weight[
i];
118 mean *= invTotWeight;
119 PRINT <<
" mean " << mean[0] <<
" " << mean[1] <<
" " << mean[2] << std::endl;
121 Eigen::Matrix3f
A = Eigen::Matrix3f::Zero();
122 for(
size_t i=0;
i<
N; ++
i) {
123 const auto diff_x = x[
i] - mean[0];
124 const auto diff_y = y[
i] - mean[1];
125 const auto diff_z = z[
i] - mean[2];
126 const auto w = weight[
i];
130 A(0,0) +=
w * diff_x * diff_x;
131 A(1,0) +=
w * diff_y * diff_x;
132 A(1,1) +=
w * diff_y * diff_y;
133 A(2,0) +=
w * diff_z * diff_x;
134 A(2,1) +=
w * diff_z * diff_y;
135 A(2,2) +=
w * diff_z * diff_z;
136 PRINT <<
" i " <<
i <<
" A " <<
A(0,0) <<
" " <<
A(0,1) <<
" " <<
A(0,2) << std::endl
137 <<
" " <<
A(1,0) <<
" " <<
A(1,1) <<
" " <<
A(1,2) << std::endl
138 <<
" " <<
A(2,0) <<
" " <<
A(2,1) <<
" " <<
A(2,2) << std::endl;
140 A *= 1./
static_cast<float>(
N);
142 PRINT <<
" A " <<
A(0,0) <<
" " <<
A(0,1) <<
" " <<
A(0,2) << std::endl
143 <<
" " <<
A(1,0) <<
" " <<
A(1,1) <<
" " <<
A(1,2) << std::endl
144 <<
" " <<
A(2,0) <<
" " <<
A(2,1) <<
" " <<
A(2,2) << std::endl;
147 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigen(A);
148 const auto& eigenValues = eigen.eigenvalues();
149 const auto& eigenVectors = eigen.eigenvectors();
152 PRINT <<
" eigenvalues " << eigenValues[0]
153 <<
" " << eigenValues[1]
154 <<
" " << eigenValues[2]
157 PRINT <<
" eigen " << eigenVectors(0,0) <<
" " << eigenVectors(0,1) <<
" " << eigenVectors(0,2) << std::endl
158 <<
" " << eigenVectors(1,0) <<
" " << eigenVectors(1,1) <<
" " << eigenVectors(1,2) << std::endl
159 <<
" " << eigenVectors(2,0) <<
" " << eigenVectors(2,1) <<
" " << eigenVectors(2,2) << std::endl;
162 auto n = eigenVectors.col(0);
163 PRINT <<
" n (unit) " <<
n[0] <<
" " <<
n[1] <<
" " <<
n[2] << std::endl;
165 const float c = -1.f * (
n[0]*mean[0] +
n[1]*mean[1] +
n[2]*mean[2]);
167 PRINT <<
" c " << c << std::endl;
170 const auto tmp = 0.5f * 1.f/
n[2];
173 const float rho2 = (1 -
sqr(
n[2]) - 4*c*
n[2]) *
sqr(tmp);
178 for(
size_t i=0;
i<
N; ++
i) {
179 const float incr =
sqr(c +
n[0]*x[
i] +
n[1]*y[
i] +
n[2]*z[
i]) * weight[
i];
180 #if defined(MK_DEBUG) || defined(EDM_ML_DEBUG) 183 <<
" chi2 incr " << incr
184 <<
" d(hit, circle) " << distance
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
double S(const TLorentzVector &, const TLorentzVector &)
std::vector< std::vector< double > > tmp
#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