CMS 3D CMS Logo

FastCircleFit.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkSeedGenerator_FastCircleFit_h
2 #define RecoTracker_TkSeedGenerator_FastCircleFit_h
3 
7 
8 #include <Eigen/Dense>
9 
10 #ifdef MK_DEBUG
11 #include <iostream>
12 #define PRINT std::cout
13 #else
15 #define PRINT LogTrace("FastCircleFit")
16 #endif
17 
18 #include <numeric>
19 #include <vector>
20 
28 public:
37  template <typename P, typename E>
38  FastCircleFit(const P& points, const E& errors) {
39  const auto N = points.size();
40  declareDynArray(float, N, x);
41  declareDynArray(float, N, y);
42  declareDynArray(float, N, z);
43  declareDynArray(float, N, weight); // 1/sigma^2
44  calculate(points, errors, x, y, z, weight);
45  }
46 
50  template <size_t N>
51  FastCircleFit(const std::array<GlobalPoint, N>& points, const std::array<GlobalError, N>& errors) {
52  std::array<float, N> x;
53  std::array<float, N> y;
54  std::array<float, N> z;
55  std::array<float, N> weight; // 1/sigma^2
56  calculate(points, errors, x, y, z, weight);
57  }
58 
59  ~FastCircleFit() = default;
60 
61  float x0() const { return x0_; }
62  float y0() const { return y0_; }
63  float rho() const { return rho_; }
64 
65  // TODO: I'm not sure if the minimized square sum is chi2 distributed
66  float chi2() const { return chi2_; }
67 
68 private:
69  template <typename P, typename E, typename C>
70  void calculate(const P& points, const E& errors, C& x, C& y, C& z, C& weight);
71 
72  template <typename T>
73  T sqr(T t) { return t*t; }
74 
75  float x0_;
76  float y0_;
77  float rho_;
78  float chi2_;
79 };
80 
81 template <typename P, typename E, typename C>
82 void FastCircleFit::calculate(const P& points, const E& errors, C& x, C& y, C& z, C& weight) {
83  const auto N = points.size();
84 
85  // transform
86  for(size_t i=0; i<N; ++i) {
87  const auto& point = points[i];
88  const auto& p = point.basicVector();
89  x[i] = p.x();
90  y[i] = p.y();
91  z[i] = sqr(p.x()) + sqr(p.y());
92 
93  // TODO: The following accounts only the hit uncertainty in phi.
94  // Albeit it "works nicely", should we account also the
95  // uncertainty in r in FPix (and the correlation between phi and r)?
96  const float phiErr2 = errors[i].phierr(point);
97  const float w = 1.f/(point.perp2()*phiErr2);
98  weight[i] = w;
99 
100  PRINT << " point " << p.x()
101  << " " << p.y()
102  << " transformed " << x[i]
103  << " " << y[i]
104  << " " << z[i]
105  << " weight " << w
106  << std::endl;
107  }
108  const float invTotWeight = 1.f/std::accumulate(weight.begin(), weight.end(), 0.f);
109  PRINT << " invTotWeight " << invTotWeight;
110 
111  Eigen::Vector3f mean = Eigen::Vector3f::Zero();
112  for(size_t i=0; i<N; ++i) {
113  const float w = weight[i];
114  mean[0] += w*x[i];
115  mean[1] += w*y[i];
116  mean[2] += w*z[i];
117  }
118  mean *= invTotWeight;
119  PRINT << " mean " << mean[0] << " " << mean[1] << " " << mean[2] << std::endl;
120 
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];
127 
128  // exploit the fact that only lower triangular part of the matrix
129  // is used in the eigenvalue calculation
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;
139  }
140  A *= 1./static_cast<float>(N);
141 
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;
145 
146  // find eigenvalues and vectors
147  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigen(A);
148  const auto& eigenValues = eigen.eigenvalues();
149  const auto& eigenVectors = eigen.eigenvectors();
150 
151  // in Eigen the first eigenvalue is the smallest
152  PRINT << " eigenvalues " << eigenValues[0]
153  << " " << eigenValues[1]
154  << " " << eigenValues[2]
155  << std::endl;
156 
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;
160 
161  // eivenvector corresponding smallest eigenvalue
162  auto n = eigenVectors.col(0);
163  PRINT << " n (unit) " << n[0] << " " << n[1] << " " << n[2] << std::endl;
164 
165  const float c = -1.f * (n[0]*mean[0] + n[1]*mean[1] + n[2]*mean[2]);
166 
167  PRINT << " c " << c << std::endl;
168 
169  // calculate fit parameters
170  const auto tmp = 0.5f * 1.f/n[2];
171  x0_ = -n[0]*tmp;
172  y0_ = -n[1]*tmp;
173  const float rho2 = (1 - sqr(n[2]) - 4*c*n[2]) * sqr(tmp);
174  rho_ = std::sqrt(rho2);
175 
176  // calculate square sum
177  float S = 0;
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)
181  const float distance = std::sqrt(sqr(x[i]-x0_) + sqr(y[i]-y0_)) - rho_;
182  PRINT << " i " << i
183  << " chi2 incr " << incr
184  << " d(hit, circle) " << distance
185  << " sigma " << 1.f/std::sqrt(weight[i])
186  << std::endl;
187 #endif
188  S += incr;
189  }
190  chi2_ = S;
191 }
192 
193 
194 #endif
float x0() const
Definition: FastCircleFit.h:61
#define PRINT
Definition: FastCircleFit.h:15
const double w
Definition: UKUtility.cc:23
~FastCircleFit()=default
FastCircleFit(const P &points, const E &errors)
Definition: FastCircleFit.h:38
Definition: weight.py:1
T sqrt(T t)
Definition: SSEVec.h:18
float y0() const
Definition: FastCircleFit.h:62
FastCircleFit(const std::array< GlobalPoint, N > &points, const std::array< GlobalError, N > &errors)
Definition: FastCircleFit.h:51
float rho() const
Definition: FastCircleFit.h:63
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
Definition: FastCircleFit.h:82
#define N
Definition: blowfish.cc:9
float chi2() const
Definition: FastCircleFit.h:66
std::pair< OmniClusterRef, TrackingParticleRef > P
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
long double T
*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
Definition: invegas.h:5