CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
FastCircleFit Class Reference

#include <FastCircleFit.h>

Public Member Functions

float chi2 () const
 
template<typename P , typename E >
 FastCircleFit (const P &points, const E &errors)
 
template<size_t N>
 FastCircleFit (const std::array< GlobalPoint, N > &points, const std::array< GlobalError, N > &errors)
 
float rho () const
 
float x0 () const
 
float y0 () const
 
 ~FastCircleFit ()=default
 

Private Member Functions

template<typename P , typename E , typename C >
void calculate (const P &points, const E &errors, C &x, C &y, C &z, C &weight)
 
template<typename T >
T sqr (T t)
 

Private Attributes

float chi2_
 
float rho_
 
float x0_
 
float y0_
 

Detailed Description

Same (almost) as FastCircle but with arbitrary number of hits.

Strandlie, Wroldsen, Frühwirth NIM A 488 (2002) 332-341. Frühwirth, Strandlie, Waltenberger NIM A 490 (2002) 366-378.

Definition at line 27 of file FastCircleFit.h.

Constructor & Destructor Documentation

template<typename P , typename E >
FastCircleFit::FastCircleFit ( const P points,
const E &  errors 
)
inline

Constructor for containers of GlobalPoint and GlobalError

Template Parameters
PContainer of GlobalPoint
EContainer of GlobalError

Container can be e.g. std::vector or DynArray.

Definition at line 38 of file FastCircleFit.h.

References calculate(), declareDynArray, N, x, y, and z.

38  {
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
45  }
Definition: weight.py:1
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
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
template<size_t N>
FastCircleFit::FastCircleFit ( const std::array< GlobalPoint, N > &  points,
const std::array< GlobalError, N > &  errors 
)
inline

Constructor for std::array of GlobalPoint and GlobalError

Definition at line 51 of file FastCircleFit.h.

References calculate(), mps_merge::weight, x, y, z, and ~FastCircleFit().

51  {
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  }
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
Definition: FastCircleFit.h:82
Definition: errors.py:1
FastCircleFit::~FastCircleFit ( )
default

Referenced by FastCircleFit().

Member Function Documentation

template<typename P , typename E , typename C >
void FastCircleFit::calculate ( const P points,
const E &  errors,
C &  x,
C &  y,
C &  z,
C &  weight 
)
private

Definition at line 82 of file FastCircleFit.h.

References patCaloMETCorrections_cff::A, EnergyCorrector::c, chi2_, SoftLeptonByDistance_cfi::distance, mps_fire::i, SiStripPI::mean, N, gen::n, AlCaHLTBitMon_ParallelJobs::p, point, PRINT, rho_, S(), sqr(), mathSSE::sqrt(), tmp, w, x0_, and y0_.

Referenced by chi2(), and FastCircleFit().

82  {
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 }
#define PRINT
Definition: FastCircleFit.h:15
const double w
Definition: UKUtility.cc:23
Definition: weight.py:1
T sqrt(T t)
Definition: SSEVec.h:18
#define N
Definition: blowfish.cc:9
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
Definition: errors.py:1
*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
float FastCircleFit::chi2 ( void  ) const
inline

Definition at line 66 of file FastCircleFit.h.

References patCaloMETCorrections_cff::C, calculate(), chi2_, hiPixelPairStep_cff::points, x, y, and z.

Referenced by CAHitQuadrupletGenerator::hitNtuplets().

66 { return chi2_; }
float FastCircleFit::rho ( ) const
inline

Definition at line 63 of file FastCircleFit.h.

References rho_.

Referenced by Lepton.Lepton::absIsoFromEA(), and Muon.Muon::absIsoWithFSR().

63 { return rho_; }
template<typename T >
T FastCircleFit::sqr ( T  t)
inlineprivate

Definition at line 73 of file FastCircleFit.h.

References lumiQTWidget::t.

Referenced by calculate().

73 { return t*t; }
float FastCircleFit::x0 ( ) const
inline

Definition at line 61 of file FastCircleFit.h.

References x0_.

61 { return x0_; }
float FastCircleFit::y0 ( ) const
inline

Definition at line 62 of file FastCircleFit.h.

References y0_.

62 { return y0_; }

Member Data Documentation

float FastCircleFit::chi2_
private

Definition at line 78 of file FastCircleFit.h.

Referenced by calculate(), and chi2().

float FastCircleFit::rho_
private

Definition at line 77 of file FastCircleFit.h.

Referenced by calculate(), and rho().

float FastCircleFit::x0_
private

Definition at line 75 of file FastCircleFit.h.

Referenced by calculate(), and x0().

float FastCircleFit::y0_
private

Definition at line 76 of file FastCircleFit.h.

Referenced by calculate(), and y0().