CMS 3D CMS Logo

Functions | Variables
WeightedMeanFitter Namespace Reference

Functions

std::pair< GlobalPoint, double > nearestPoint (const GlobalPoint &vertex, reco::Track iclus)
 
TransientVertex weightedMeanOutlierRejection (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus)
 
TransientVertex weightedMeanOutlierRejectionBeamSpot (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus, const reco::BeamSpot &beamSpot)
 
TransientVertex weightedMeanOutlierRejectionVarianceAsError (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< std::vector< reco::TransientTrack >>::const_iterator iclus)
 

Variables

constexpr float corr_x = 1.2
 
constexpr float corr_x_bs = 1.0
 
constexpr float corr_z = 1.4
 
constexpr int maxIterations = 50
 
constexpr float muSquare = 9.
 
constexpr float precision = 1e-24
 
constexpr float startError = 20.0
 

Function Documentation

◆ nearestPoint()

std::pair<GlobalPoint, double> WeightedMeanFitter::nearestPoint ( const GlobalPoint vertex,
reco::Track  iclus 
)
inline

Definition at line 21 of file WeightedMeanFitter.h.

References AlCaHLTBitMon_ParallelJobs::p, funct::pow(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), mathSSE::sqrt(), submitPVValidationJobs::t, bphysicsOniaDQM_cfi::vertex, btvMC_cff::vx, reco::TrackBase::vx(), btvMC_cff::vy, reco::TrackBase::vy(), btvMC_cff::vz, and reco::TrackBase::vz().

Referenced by weightedMeanOutlierRejection(), and weightedMeanOutlierRejectionBeamSpot().

21  {
22  double ox = iclus.vx();
23  double oy = iclus.vy();
24  double oz = iclus.vz();
25 
26  double vx = iclus.px();
27  double vy = iclus.py();
28  double vz = iclus.pz();
29 
30  double opx = vertex.x() - ox;
31  double opy = vertex.y() - oy;
32  double opz = vertex.z() - oz;
33 
34  double vnorm2 = (vx * vx + vy * vy + vz * vz);
35  double t = (vx * opx + vy * opy + vz * opz) / (vnorm2);
36 
37  GlobalPoint p(ox + t * vx, oy + t * vy, oz + t * vz);
38  return std::pair<GlobalPoint, double>(
39  p,
40  std::sqrt(std::pow(p.x() - vertex.x(), 2) + std::pow(p.y() - vertex.y(), 2) + std::pow(p.z() - vertex.z(), 2)));
41  }
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:655
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:640
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:643
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
T sqrt(T t)
Definition: SSEVec.h:19
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ weightedMeanOutlierRejection()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejection ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< reco::TransientTrack iclus 
)
inline

Definition at line 43 of file WeightedMeanFitter.h.

References funct::abs(), nano_mu_local_reco_cff::chi2, corr_x, corr_z, submitPVResolutionJobs::err, err_x, err_z, mps_fire::i, createfilelist::int, muSquare, ndof, nearestPoint(), AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, funct::pow(), precision, edm::second(), startError, HLT_2024v13_cff::track, findQualityFiles::v, and x.

Referenced by WeightedMeanPrimaryVertexEstimator::fit().

44  {
45  float x = 0., y = 0., z = 0.;
46  float s_wx = 0., s_wz = 0.;
47  float s2_wx = 0., s2_wz = 0.;
48  float wx = 0., wz = 0., chi2 = 0.;
49  float ndof_x = 0.;
50 
52  err(0, 0) = startError / 10 * startError / 10;
53  err(1, 1) = startError / 10 * startError / 10;
54  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
55  for (const auto& p : points) {
56  wx = p.second.x() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.x(), 2);
57 
58  wz = p.second.z() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.z(), 2);
59 
60  x += p.first.x() * wx;
61  y += p.first.y() * wx;
62  z += p.first.z() * wz;
63 
64  s_wx += wx;
65  s_wz += wz;
66  }
67 
68  if (s_wx == 0. || s_wz == 0.) {
69  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
70  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
71  }
72 
73  x /= s_wx;
74  y /= s_wx;
75  z /= s_wz;
76 
77  float old_x, old_y, old_z;
78 
79  float xpull;
80  int niter = 0;
81 
82  float err_x, err_z;
83 
84  err_x = 1. / s_wx;
85  err_z = 1. / s_wz;
86 
87  while ((niter++) < 2) {
88  old_x = x;
89  old_y = y;
90  old_z = z;
91  s_wx = 0;
92  s_wz = 0;
93  s2_wx = 0;
94  s2_wz = 0;
95  x = 0;
96  y = 0;
97  z = 0;
98  ndof_x = 0;
99 
100  for (unsigned int i = 0; i < (unsigned int)points.size(); i++) {
101  std::pair<GlobalPoint, double> p = nearestPoint(GlobalPoint(old_x, old_y, old_z), (iclus)[i].track());
102 
103  wx = points[i].second.x() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.x(), 2);
104  wz = points[i].second.z() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.z(), 2);
105 
106  xpull = 0.;
107 
108  if (std::pow(p.first.x() - old_x, 2) / (wx + err_x) < muSquare &&
109  std::pow(p.first.y() - old_y, 2) / (wx + err_x) < muSquare &&
110  std::pow(p.first.z() - old_z, 2) / (wz + err_z) < muSquare)
111  xpull = 1.;
112 
113  ndof_x += xpull;
114 
115  wx = xpull / wx;
116  wz = xpull / wz;
117 
118  x += wx * p.first.x();
119  y += wx * p.first.y();
120  z += wz * p.first.z();
121 
122  s_wx += wx;
123  s_wz += wz;
124 
125  s2_wx += wx * xpull;
126  s2_wz += wz * xpull;
127  }
128 
129  if (s_wx == 0. || s_wz == 0.) {
130  edm::LogWarning("WeightedMeanFitter")
131  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
132  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
133  }
134  x /= s_wx;
135  y /= s_wx;
136  z /= s_wz;
137 
138  err_x = (s2_wx / std::pow(s_wx, 2));
139  err_z = (s2_wz / std::pow(s_wz, 2));
140 
141  if (std::abs(x - old_x) < (precision / 1.) && std::abs(y - old_y) < (precision / 1.) &&
142  std::abs(z - old_z) < (precision / 1.)) {
143  break;
144  }
145  }
146 
147  err(0, 0) = err_x * corr_x * corr_x;
148  err(1, 1) = err_x * corr_x * corr_x;
149  err(2, 2) = err_z * corr_z * corr_z;
150 
151  float dist = 0;
152  for (const auto& p : points) {
153  wx = p.second.x();
154  wx = wx <= precision ? precision : wx;
155 
156  wz = p.second.z();
157  wz = wz <= precision ? precision : wz;
158 
159  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + err(0, 0));
160  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + err(1, 1));
161  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + err(2, 2));
162  chi2 += dist;
163  }
164  float ndof =
165  ndof_x > 1 ? (2 * ndof_x - 3) : 0.00001; // ndof_x is actually the number of tracks with non-zero weight
166  TransientVertex v(GlobalPoint(x, y, z), err, iclus, chi2, ndof);
167  return v;
168  }
constexpr float corr_x
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
float float float z
U second(std::pair< T, U > const &p)
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr float corr_z
std::pair< GlobalPoint, double > nearestPoint(const GlobalPoint &vertex, reco::Track iclus)
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Log< level::Warning, false > LogWarning
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
constexpr float startError

◆ weightedMeanOutlierRejectionBeamSpot()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejectionBeamSpot ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< reco::TransientTrack iclus,
const reco::BeamSpot beamSpot 
)
inline

Definition at line 170 of file WeightedMeanFitter.h.

References funct::abs(), pwdgSkimBPark_cfi::beamSpot, nano_mu_local_reco_cff::chi2, corr_x_bs, corr_z, submitPVResolutionJobs::err, err_x, err_y, err_z, mps_fire::i, createfilelist::int, muSquare, nearestPoint(), AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, funct::pow(), precision, edm::second(), startError, HLT_2024v13_cff::track, findQualityFiles::v, and x.

Referenced by WeightedMeanPrimaryVertexEstimator::fit().

173  {
174  float x = 0., y = 0., z = 0.;
175  float s_wx = 0., s_wz = 0.;
176  float s2_wx = 0., s2_wz = 0.;
177  float wx = 0., wz = 0., chi2 = 0.;
178  float wy = 0., s_wy = 0., s2_wy = 0.;
179  float ndof_x = 0.;
180 
182  err(0, 0) = startError / 10 * startError / 10;
183  err(1, 1) = startError / 10 * startError / 10;
184  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
185 
186  GlobalError bse(beamSpot.rotatedCovariance3D());
187  GlobalPoint bsp(Basic3DVector<float>(beamSpot.position()));
188 
189  for (const auto& p : points) {
190  wx = p.second.x() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.x(), 2);
191  wy = p.second.y() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.y(), 2);
192 
193  wz = p.second.z() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.z(), 2);
194 
195  x += p.first.x() * wx;
196  y += p.first.y() * wy;
197  z += p.first.z() * wz;
198 
199  s_wx += wx;
200  s_wy += wy;
201  s_wz += wz;
202  }
203 
204  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
205  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
206  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
207  }
208  // use the square of covariance element to increase it's weight: it will be the most important
209  wx = bse.cxx() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(bse.cxx(), 2);
210  wy = bse.cyy() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(bse.cyy(), 2);
211 
212  x += bsp.x() * wx;
213  y += bsp.y() * wy;
214 
215  x /= (s_wx + wx);
216  y /= (s_wy + wy);
217  z /= s_wz;
218 
219  float old_x, old_y, old_z;
220 
221  float xpull;
222  int niter = 0;
223 
224  float err_x, err_y, err_z;
225 
226  err_x = 1. / s_wx;
227  err_y = 1. / s_wy;
228  err_z = 1. / s_wz;
229 
230  while ((niter++) < 2) {
231  old_x = x;
232  old_y = y;
233  old_z = z;
234  s_wx = 0;
235  s_wz = 0;
236  s2_wx = 0;
237  s2_wz = 0;
238 
239  s_wy = 0;
240  s2_wy = 0;
241 
242  x = 0;
243  y = 0;
244  z = 0;
245  ndof_x = 0;
246 
247  for (unsigned int i = 0; i < (unsigned int)points.size(); i++) {
248  std::pair<GlobalPoint, double> p = nearestPoint(GlobalPoint(old_x, old_y, old_z), (iclus)[i].track());
249 
250  wx = points[i].second.x() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.x(), 2);
251  wy = points[i].second.y() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.y(), 2);
252 
253  wz = points[i].second.z() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.z(), 2);
254 
255  xpull = 0.;
256  if (std::pow(p.first.x() - old_x, 2) / (wx + err_x) < muSquare &&
257  std::pow(p.first.y() - old_y, 2) / (wy + err_y) < muSquare &&
258  std::pow(p.first.z() - old_z, 2) / (wz + err_z) < muSquare)
259  xpull = 1.;
260 
261  ndof_x += xpull;
262 
263  wx = xpull / wx;
264  wy = xpull / wy;
265  wz = xpull / wz;
266 
267  x += wx * p.first.x();
268  y += wy * p.first.y();
269  z += wz * p.first.z();
270 
271  s_wx += wx;
272  s_wy += wy;
273  s_wz += wz;
274 
275  s2_wx += wx * xpull;
276  s2_wy += wy * xpull;
277  s2_wz += wz * xpull;
278  }
279 
280  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
281  edm::LogWarning("WeightedMeanFitter")
282  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
283  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
284  }
285  wx = bse.cxx() <= std::pow(precision, 2) ? 1. / std::pow(precision, 2) : 1. / bse.cxx();
286  wy = bse.cyy() <= std::pow(precision, 2) ? 1. / std::pow(precision, 2) : 1. / bse.cyy();
287 
288  x += bsp.x() * wx;
289  y += bsp.y() * wy;
290  s_wx += wx;
291  s2_wx += wx;
292  s_wy += wy;
293  s2_wy += wy;
294 
295  x /= s_wx;
296  y /= s_wy;
297  z /= s_wz;
298 
299  err_x = (s2_wx / std::pow(s_wx, 2));
300  err_y = (s2_wy / std::pow(s_wy, 2));
301  err_z = (s2_wz / std::pow(s_wz, 2));
302 
303  if (std::abs(x - old_x) < (precision) && std::abs(y - old_y) < (precision) && std::abs(z - old_z) < (precision)) {
304  break;
305  }
306  }
307  err(0, 0) = err_x * corr_x_bs * corr_x_bs;
308  err(1, 1) = err_y * corr_x_bs * corr_x_bs;
309  err(2, 2) = err_z * corr_z * corr_z;
310 
311  float dist = 0;
312  for (const auto& p : points) {
313  wx = p.second.x();
314  wx = wx <= precision ? precision : wx;
315 
316  wz = p.second.z();
317  wz = wz <= precision ? precision : wz;
318 
319  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + err(0, 0));
320  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + err(1, 1));
321  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + err(2, 2));
322  chi2 += dist;
323  }
324  TransientVertex v(GlobalPoint(x, y, z), err, iclus, chi2, (int)ndof_x);
325  return v;
326  }
constexpr float corr_x_bs
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
float float float z
U second(std::pair< T, U > const &p)
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr float corr_z
std::pair< GlobalPoint, double > nearestPoint(const GlobalPoint &vertex, reco::Track iclus)
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
float err_y
Log< level::Warning, false > LogWarning
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
constexpr float startError

◆ weightedMeanOutlierRejectionVarianceAsError()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejectionVarianceAsError ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< std::vector< reco::TransientTrack >>::const_iterator  iclus 
)
inline

Definition at line 328 of file WeightedMeanFitter.h.

References funct::abs(), nano_mu_local_reco_cff::chi2, submitPVResolutionJobs::err, err_x, err_y, err_z, JetChargeProducer_cfi::exp, maxIterations, muSquare, AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, funct::pow(), precision, startError, findQualityFiles::v, and x.

330  {
331  float x = 0, y = 0, z = 0, s_wx = 0, s_wy = 0, s_wz = 0, wx = 0, wy = 0, wz = 0, chi2 = 0;
332  float ndof_x = 0;
334  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
335  err(0, 0) = err(1, 1) = err(2, 2) / 100.;
336 
337  for (const auto& p : points) {
338  wx = p.second.x();
339  wx = wx <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(wx, 2);
340 
341  wz = p.second.z();
342  wz = wz <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(wz, 2);
343 
344  x += p.first.x() * wx;
345  y += p.first.y() * wx;
346  z += p.first.z() * wz;
347 
348  s_wx += wx;
349  s_wz += wz;
350  }
351 
352  if (s_wx == 0. || s_wz == 0.) {
353  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
354  return TransientVertex(GlobalPoint(0, 0, 0), err, *iclus, 0, 0);
355  }
356 
357  x /= s_wx;
358  y /= s_wx;
359  z /= s_wz;
360 
361  float old_x, old_y, old_z;
362  float xpull;
363  int niter = 0;
364  float err_x, err_y, err_z;
365  err_x = 1. / s_wx;
366  err_y = 1. / s_wx;
367  err_z = 1. / s_wz;
368  float s_err_x = 0., s_err_y = 0., s_err_z = 0.;
369  while ((niter++) < maxIterations) {
370  old_x = x;
371  old_y = y;
372  old_z = z;
373  s_wx = 0;
374  s_wy = 0;
375  s_wz = 0;
376  x = 0;
377  y = 0;
378  z = 0;
379  s_err_x = 0.;
380  s_err_y = 0.;
381  s_err_z = 0.;
382 
383  for (const auto& p : points) {
384  wx = p.second.x();
385  wx = wx <= precision ? precision : wx;
386 
387  wy = wx * wx + err_y;
388  wx = wx * wx + err_x;
389 
390  wz = p.second.z();
391  wz = wz <= precision ? precision : wz;
392  wz = wz * wz + err_z;
393 
394  xpull = std::pow((p.first.x() - old_x), 2) / wx;
395  xpull += std::pow((p.first.y() - old_y), 2) / wy;
396  xpull += std::pow((p.first.z() - old_z), 2) / wz;
397  xpull = 1. / (1. + std::exp(-0.5 * ((muSquare)-xpull)));
398  ndof_x += xpull;
399 
400  wx = 1. / wx;
401  wy = 1. / wy;
402  wz = 1. / wz;
403 
404  wx *= xpull;
405  wy *= xpull;
406  wz *= xpull;
407 
408  x += wx * p.first.x();
409  y += wy * p.first.y();
410  z += wz * p.first.z();
411 
412  s_wx += wx;
413  s_wy += wy;
414  s_wz += wz;
415 
416  s_err_x += wx * pow(p.first.x() - old_x, 2);
417  s_err_y += wy * pow(p.first.y() - old_y, 2);
418  s_err_z += wz * pow(p.first.z() - old_z, 2);
419  }
420  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
421  edm::LogWarning("WeightedMeanFitter")
422  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
423  return TransientVertex(GlobalPoint(0, 0, 0), err, *iclus, 0, 0);
424  }
425  x /= s_wx;
426  y /= s_wy;
427  z /= s_wz;
428 
429  err_x = s_err_x / s_wx;
430  err_y = s_err_y / s_wy;
431  err_z = s_err_z / s_wz;
432 
433  if (std::abs(x - old_x) < (precision / 1.) && std::abs(y - old_y) < (precision / 1.) &&
434  std::abs(z - old_z) < (precision / 1.))
435  break;
436  }
437 
438  err(0, 0) = err_x;
439  err(1, 1) = err_y;
440  err(2, 2) = err_z;
441 
442  float dist = 0.f;
443  for (const auto& p : points) {
444  wx = p.second.x();
445  wx = wx <= precision ? precision : wx;
446 
447  wz = p.second.z();
448  wz = wz <= precision ? precision : wz;
449 
450  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + std::pow(err(0, 0), 2));
451  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + std::pow(err(1, 1), 2));
452  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + std::pow(err(2, 2), 2));
453  chi2 += dist;
454  }
455  TransientVertex v(GlobalPoint(x, y, z), err, *iclus, chi2, (int)ndof_x);
456  return v;
457  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
float float float z
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr int maxIterations
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
float err_y
Log< level::Warning, false > LogWarning
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
constexpr float startError

Variable Documentation

◆ corr_x

constexpr float WeightedMeanFitter::corr_x = 1.2

Definition at line 15 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejection().

◆ corr_x_bs

constexpr float WeightedMeanFitter::corr_x_bs = 1.0

Definition at line 16 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejectionBeamSpot().

◆ corr_z

constexpr float WeightedMeanFitter::corr_z = 1.4

◆ maxIterations

constexpr int WeightedMeanFitter::maxIterations = 50

Definition at line 18 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejectionVarianceAsError().

◆ muSquare

constexpr float WeightedMeanFitter::muSquare = 9.

◆ precision

constexpr float WeightedMeanFitter::precision = 1e-24

◆ startError

constexpr float WeightedMeanFitter::startError = 20.0