CMS 3D CMS Logo

MuonResidualsFitter.cc
Go to the documentation of this file.
1 // $Id: MuonResidualsFitter.cc,v 1.11 2011/04/15 21:48:21 khotilov Exp $
2 
3 #ifdef STANDALONE_FITTER
4 #include "MuonResidualsFitter.h"
5 #else
7 #endif
8 
9 #include <fstream>
10 #include <set>
11 #include "TMath.h"
12 #include "TH1.h"
13 #include "TF1.h"
14 #include "TRobustEstimator.h"
15 
16 // all global variables begin with "MuonResidualsFitter_" to avoid
17 // namespace clashes (that is, they do what would ordinarily be done
18 // with a class structure, but Minuit requires them to be global)
19 
20 const double MuonResidualsFitter_gsbinsize = 0.01;
21 const double MuonResidualsFitter_tsbinsize = 0.1;
24 
27 
29 
30 // fit function
31 double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma) {
32  sigma = fabs(sigma);
33  static const double cgaus = 0.5 * log(2. * M_PI);
34  return (-pow(residual - center, 2) * 0.5 / sigma / sigma) - cgaus - log(sigma);
35 }
36 
37 // TF1 interface version
38 Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par) {
39  return par[0] * exp(MuonResidualsFitter_logPureGaussian(xvec[0], par[1], par[2]));
40 }
41 
42 // 2D Gauss fit function
43 double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r) {
44  sx = fabs(sx);
45  sy = fabs(sy);
46  static const double c2gaus = log(2. * M_PI);
47  double one_r2 = 1. - r * r;
48  double dx = x - x0;
49  double dy = y - y0;
50  return -0.5 / one_r2 * (pow(dx / sx, 2) + pow(dy / sy, 2) - 2. * r * dx / sx * dy / sy) - c2gaus -
51  log(sx * sy * sqrt(one_r2));
52 }
53 
55  double toversigma, double gammaoversigma, double max, double step, double power) {
56  if (gammaoversigma == 0.)
57  return (-toversigma * toversigma / 2.) - log(sqrt(2 * M_PI));
58 
59  double sum = 0.;
60  double uplus = 0.;
61  double integrandplus_last = 0.;
62  double integrandminus_last = 0.;
63  for (double inc = 0.; uplus < max; inc += step) {
64  double uplus_last = uplus;
65  uplus = pow(inc, power);
66 
67  double integrandplus = exp(-pow(uplus - toversigma, 2) / 2.) / (uplus * uplus / gammaoversigma + gammaoversigma) /
68  M_PI / sqrt(2. * M_PI);
69  double integrandminus = exp(-pow(-uplus - toversigma, 2) / 2.) / (uplus * uplus / gammaoversigma + gammaoversigma) /
70  M_PI / sqrt(2. * M_PI);
71 
72  sum += integrandplus * (uplus - uplus_last);
73  sum += 0.5 * fabs(integrandplus_last - integrandplus) * (uplus - uplus_last);
74 
75  sum += integrandminus * (uplus - uplus_last);
76  sum += 0.5 * fabs(integrandminus_last - integrandminus) * (uplus - uplus_last);
77 
78  integrandplus_last = integrandplus;
79  integrandminus_last = integrandminus;
80  }
81  return log(sum);
82 }
83 
84 // fit function
85 double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma) {
86  sigma = fabs(sigma);
87  double gammaoversigma = fabs(gamma / sigma);
88  double toversigma = fabs((residual - center) / sigma);
89 
90  int gsbin0 = int(floor(gammaoversigma / MuonResidualsFitter_gsbinsize));
91  int gsbin1 = int(ceil(gammaoversigma / MuonResidualsFitter_gsbinsize));
92  int tsbin0 = int(floor(toversigma / MuonResidualsFitter_tsbinsize));
93  int tsbin1 = int(ceil(toversigma / MuonResidualsFitter_tsbinsize));
94 
95  bool gsisbad = (gsbin0 >= MuonResidualsFitter_numgsbins || gsbin1 >= MuonResidualsFitter_numgsbins);
96  bool tsisbad = (tsbin0 >= MuonResidualsFitter_numtsbins || tsbin1 >= MuonResidualsFitter_numtsbins);
97 
98  if (gsisbad || tsisbad) {
99  return log(gammaoversigma / M_PI) - log(toversigma * toversigma + gammaoversigma * gammaoversigma) - log(sigma);
100  } else {
101  double val00 = MuonResidualsFitter_lookup_table[gsbin0][tsbin0];
102  double val01 = MuonResidualsFitter_lookup_table[gsbin0][tsbin1];
103  double val10 = MuonResidualsFitter_lookup_table[gsbin1][tsbin0];
104  double val11 = MuonResidualsFitter_lookup_table[gsbin1][tsbin1];
105 
106  double val0 = val00 + ((toversigma / MuonResidualsFitter_tsbinsize) - tsbin0) * (val01 - val00);
107  double val1 = val10 + ((toversigma / MuonResidualsFitter_tsbinsize) - tsbin0) * (val11 - val10);
108 
109  double val = val0 + ((gammaoversigma / MuonResidualsFitter_gsbinsize) - gsbin0) * (val1 - val0);
110 
111  return val - log(sigma);
112  }
113 }
114 
115 // TF1 interface version
116 Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par) {
117  return par[0] * exp(MuonResidualsFitter_logPowerLawTails(xvec[0], par[1], par[2], par[3]));
118 }
119 
120 double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma) {
121  return log(TMath::Voigt(residual - center, fabs(sigma), fabs(gamma) * 2.));
122 }
123 
124 Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par) {
125  return par[0] * TMath::Voigt(xvec[0] - par[1], fabs(par[2]), fabs(par[3]) * 2.);
126 }
127 
128 double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma) {
129  double x = residual - center;
130  double s = fabs(sigma);
131  double m = 2. * s;
132  static const double e2 = exp(-2.), sqerf = sqrt(2. * M_PI) * erf(sqrt(2.));
133  double a = pow(m, 4) * e2;
134  double n = sqerf * s + 2. * a * pow(m, -3) / 3.;
135 
136  if (fabs(x) < m)
137  return -x * x / (2. * s * s) - log(n);
138  else
139  return log(a) - 4. * log(fabs(x)) - log(n);
140 }
141 
142 Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par) {
143  return par[0] * exp(MuonResidualsFitter_logGaussPowerTails(xvec[0], par[1], par[2]));
144 }
145 
146 double MuonResidualsFitter_integrate_pureGaussian(double low, double high, double center, double sigma) {
147  static const double isqr2 = 1. / sqrt(2.);
148  return (erf((high + center) * isqr2 / sigma) - erf((low + center) * isqr2 / sigma)) * exp(0.5 / sigma / sigma) * 0.5;
149 }
150 
152  : m_residualsModel(residualsModel),
153  m_minHits(minHits),
154  m_useResiduals(useResiduals),
155  m_weightAlignment(weightAlignment),
156  m_printLevel(0),
157  m_strategy(1),
158  m_cov(1),
159  m_loglikelihood(0.) {
162  throw cms::Exception("MuonResidualsFitter") << "unrecognized residualsModel";
163 }
164 
166  for (std::vector<double *>::const_iterator residual = residuals_begin(); residual != residuals_end(); ++residual) {
167  delete[](*residual);
168  }
169 }
170 
171 void MuonResidualsFitter::fix(int parNum, bool dofix) {
172  assert(0 <= parNum && parNum < npar());
173  if (m_fixed.empty())
174  m_fixed.resize(npar(), false);
175  m_fixed[parNum] = dofix;
176 }
177 
178 bool MuonResidualsFitter::fixed(int parNum) {
179  assert(0 <= parNum && parNum < npar());
180  if (m_fixed.empty())
181  return false;
182  else
183  return m_fixed[parNum];
184 }
185 
186 void MuonResidualsFitter::fill(double *residual) {
187  m_residuals.push_back(residual);
188  m_residuals_ok.push_back(true);
189 }
190 
191 double MuonResidualsFitter::covarianceElement(int parNum1, int parNum2) {
192  assert(0 <= parNum1 && parNum1 < npar());
193  assert(0 <= parNum2 && parNum2 < npar());
194  assert(m_cov.GetNcols() == npar()); // m_cov might have not yet been resized to account for proper #parameters
195  return m_cov(parNum2parIdx(parNum1), parNum2parIdx(parNum2));
196 }
197 
199  long num = 0;
200  for (std::vector<double *>::const_iterator resiter = residuals_begin(); resiter != residuals_end(); ++resiter)
201  num++;
202  return num;
203 }
204 
207  return;
209 
210  std::ifstream convolution_table("convolution_table.txt");
211  if (convolution_table.is_open()) {
212  int numgsbins = 0;
213  int numtsbins = 0;
214  double tsbinsize = 0.;
215  double gsbinsize = 0.;
216 
217  convolution_table >> numgsbins >> numtsbins >> tsbinsize >> gsbinsize;
218  if (numgsbins != MuonResidualsFitter_numgsbins || numtsbins != MuonResidualsFitter_numtsbins ||
219  tsbinsize != MuonResidualsFitter_tsbinsize || gsbinsize != MuonResidualsFitter_gsbinsize) {
220  throw cms::Exception("MuonResidualsFitter") << "convolution_table.txt has the wrong bin width/bin size. Throw "
221  "it away and let the fitter re-create the file.\n";
222  }
223 
224  for (int gsbin = 0; gsbin < MuonResidualsFitter_numgsbins; gsbin++) {
225  for (int tsbin = 0; tsbin < MuonResidualsFitter_numtsbins; tsbin++) {
226  int read_gsbin = 0;
227  int read_tsbin = 0;
228  double val = 0.;
229 
230  convolution_table >> read_gsbin >> read_tsbin >> val;
231  if (read_gsbin != gsbin || read_tsbin != tsbin) {
232  throw cms::Exception("MuonResidualsFitter")
233  << "convolution_table.txt is out of order. Throw it away and let the fitter re-create the file.\n";
234  }
235  MuonResidualsFitter_lookup_table[gsbin][tsbin] = val;
236  }
237  }
238  convolution_table.close();
239  } else {
240  std::ofstream convolution_table2("convolution_table.txt");
241 
242  if (!convolution_table2.is_open())
243  throw cms::Exception("MuonResidualsFitter") << "Couldn't write to file convolution_table.txt\n";
244 
245  convolution_table2 << MuonResidualsFitter_numgsbins << " " << MuonResidualsFitter_numtsbins << " "
247 
248  std::cout << "Initializing convolution look-up table (takes a few minutes)..." << std::endl;
249 
250  for (int gsbin = 0; gsbin < MuonResidualsFitter_numgsbins; gsbin++) {
251  double gammaoversigma = double(gsbin) * MuonResidualsFitter_gsbinsize;
252  std::cout << " gsbin " << gsbin << "/" << MuonResidualsFitter_numgsbins << std::endl;
253  for (int tsbin = 0; tsbin < MuonResidualsFitter_numtsbins; tsbin++) {
254  double toversigma = double(tsbin) * MuonResidualsFitter_tsbinsize;
255 
256  // 1e-6 errors (out of a value of ~0.01) with max=100, step=0.001, power=4 (max=1000 does a little better with the tails)
257  MuonResidualsFitter_lookup_table[gsbin][tsbin] =
258  MuonResidualsFitter_compute_log_convolution(toversigma, gammaoversigma);
259 
260  // <10% errors with max=20, step=0.005, power=4 (faster computation for testing)
261  // MuonResidualsFitter_lookup_table[gsbin][tsbin] = MuonResidualsFitter_compute_log_convolution(toversigma, gammaoversigma, 100., 0.005, 4.);
262 
263  convolution_table2 << gsbin << " " << tsbin << " " << MuonResidualsFitter_lookup_table[gsbin][tsbin]
264  << std::endl;
265  }
266  }
267  convolution_table2.close();
268  std::cout << "Initialization done!" << std::endl;
269  }
270 }
271 
272 bool MuonResidualsFitter::dofit(void (*fcn)(int &, double *, double &, double *, int),
273  std::vector<int> &parNum,
274  std::vector<std::string> &parName,
275  std::vector<double> &start,
276  std::vector<double> &step,
277  std::vector<double> &low,
278  std::vector<double> &high) {
280 
281  MuonResidualsFitter_TMinuit = new TMinuit(npar());
282  // MuonResidualsFitter_TMinuit->SetPrintLevel(m_printLevel);
283  MuonResidualsFitter_TMinuit->SetPrintLevel();
284  MuonResidualsFitter_TMinuit->SetObjectFit(fitinfo);
287 
288  std::vector<int>::const_iterator iNum = parNum.begin();
289  std::vector<std::string>::const_iterator iName = parName.begin();
290  std::vector<double>::const_iterator istart = start.begin();
291  std::vector<double>::const_iterator istep = step.begin();
292  std::vector<double>::const_iterator ilow = low.begin();
293  std::vector<double>::const_iterator ihigh = high.begin();
294 
295  //MuonResidualsFitter_TMinuit->SetPrintLevel(-1);
296 
297  for (; iNum != parNum.end(); ++iNum, ++iName, ++istart, ++istep, ++ilow, ++ihigh) {
298  MuonResidualsFitter_TMinuit->DefineParameter(*iNum, iName->c_str(), *istart, *istep, *ilow, *ihigh);
299  if (fixed(*iNum))
300  MuonResidualsFitter_TMinuit->FixParameter(*iNum);
301  }
302 
303  double arglist[10];
304  int ierflg;
305  int smierflg; //second MIGRAD ierflg
306 
307  // chi^2 errors should be 1.0, log-likelihood should be 0.5
308  for (int i = 0; i < 10; i++)
309  arglist[i] = 0.;
310  arglist[0] = 0.5;
311  ierflg = 0;
312  smierflg = 0;
313  MuonResidualsFitter_TMinuit->mnexcm("SET ERR", arglist, 1, ierflg);
314  if (ierflg != 0) {
316  delete fitinfo;
317  return false;
318  }
319 
320  // set strategy = 2 (more refined fits)
321  for (int i = 0; i < 10; i++)
322  arglist[i] = 0.;
323  arglist[0] = m_strategy;
324  ierflg = 0;
325  MuonResidualsFitter_TMinuit->mnexcm("SET STR", arglist, 1, ierflg);
326  if (ierflg != 0) {
328  delete fitinfo;
329  return false;
330  }
331 
332  bool try_again = false;
333 
334  // minimize
335  for (int i = 0; i < 10; i++)
336  arglist[i] = 0.;
337  arglist[0] = 50000;
338  ierflg = 0;
339  MuonResidualsFitter_TMinuit->mnexcm("MIGRAD", arglist, 1, ierflg);
340  if (ierflg != 0)
341  try_again = true;
342 
343  // just once more, if needed (using the final Minuit parameters from the failed fit; often works)
344  if (try_again) {
345  for (int i = 0; i < 10; i++)
346  arglist[i] = 0.;
347  arglist[0] = 50000;
348  MuonResidualsFitter_TMinuit->mnexcm("MIGRAD", arglist, 1, smierflg);
349  }
350 
351  Double_t fmin, fedm, errdef;
352  Int_t npari, nparx, istat;
353  MuonResidualsFitter_TMinuit->mnstat(fmin, fedm, errdef, npari, nparx, istat);
354 
355  if (istat != 3) {
356  for (int i = 0; i < 10; i++)
357  arglist[i] = 0.;
358  ierflg = 0;
359  MuonResidualsFitter_TMinuit->mnexcm("HESSE", arglist, 0, ierflg);
360  }
361 
362  // read-out the results
363  m_loglikelihood = -fmin;
364 
365  m_value.clear();
366  m_error.clear();
367  for (int i = 0; i < npar(); i++) {
368  double v, e;
369  MuonResidualsFitter_TMinuit->GetParameter(i, v, e);
370  m_value.push_back(v);
371  m_error.push_back(e);
372  }
373  m_cov.ResizeTo(npar(), npar());
374  MuonResidualsFitter_TMinuit->mnemat(m_cov.GetMatrixArray(), npar());
375 
377  delete fitinfo;
378  if (smierflg != 0)
379  return false;
380  return true;
381 }
382 
384  long rows = numResiduals();
385  int cols = ndata();
386  int whichcopy = which;
387 
388  fwrite(&rows, sizeof(long), 1, file);
389  fwrite(&cols, sizeof(int), 1, file);
390  fwrite(&whichcopy, sizeof(int), 1, file);
391 
392  double *likeAChecksum = new double[cols];
393  double *likeAChecksum2 = new double[cols];
394  for (int i = 0; i < cols; i++) {
395  likeAChecksum[i] = 0.;
396  likeAChecksum2[i] = 0.;
397  }
398 
399  for (std::vector<double *>::const_iterator residual = residuals_begin(); residual != residuals_end(); ++residual) {
400  fwrite((*residual), sizeof(double), cols, file);
401  for (int i = 0; i < cols; i++) {
402  if (fabs((*residual)[i]) > likeAChecksum[i])
403  likeAChecksum[i] = fabs((*residual)[i]);
404  if (fabs((*residual)[i]) < likeAChecksum2[i])
405  likeAChecksum2[i] = fabs((*residual)[i]);
406  }
407  } // end loop over residuals
408 
409  // the idea is that mal-formed doubles are likely to be huge values (or tiny values)
410  // because the exponent gets screwed up; we want to check for that
411  fwrite(likeAChecksum, sizeof(double), cols, file);
412  fwrite(likeAChecksum2, sizeof(double), cols, file);
413 
414  delete[] likeAChecksum;
415  delete[] likeAChecksum2;
416 }
417 
419  long rows = -100;
420  int cols = -100;
421  int readwhich = -100;
422 
423  fread(&rows, sizeof(long), 1, file);
424  fread(&cols, sizeof(int), 1, file);
425  fread(&readwhich, sizeof(int), 1, file);
426 
427  if (cols != ndata() || rows < 0 || readwhich != which) {
428  throw cms::Exception("MuonResidualsFitter")
429  << "temporary file is corrupted (which = " << which << " readwhich = " << readwhich << " rows = " << rows
430  << " cols = " << cols << ")\n";
431  }
432 
433  double *likeAChecksum = new double[cols];
434  double *likeAChecksum2 = new double[cols];
435  for (int i = 0; i < cols; i++) {
436  likeAChecksum[i] = 0.;
437  likeAChecksum2[i] = 0.;
438  }
439 
440  m_residuals.reserve(rows);
441  for (long row = 0; row < rows; row++) {
442  double *residual = new double[cols];
443  fread(residual, sizeof(double), cols, file);
444  fill(residual);
445  for (int i = 0; i < cols; i++) {
446  if (fabs(residual[i]) > likeAChecksum[i])
447  likeAChecksum[i] = fabs(residual[i]);
448  if (fabs(residual[i]) < likeAChecksum2[i])
449  likeAChecksum2[i] = fabs(residual[i]);
450  }
451  } // end loop over records in file
452 
453  double *readChecksum = new double[cols];
454  double *readChecksum2 = new double[cols];
455  fread(readChecksum, sizeof(double), cols, file);
456  fread(readChecksum2, sizeof(double), cols, file);
457 
458  for (int i = 0; i < cols; i++) {
459  if (fabs(likeAChecksum[i] - readChecksum[i]) > 1e-10 ||
460  fabs(1. / likeAChecksum2[i] - 1. / readChecksum2[i]) > 1e10) {
461  throw cms::Exception("MuonResidualsFitter")
462  << "temporary file is corrupted (which = " << which << " rows = " << rows << " likeAChecksum "
463  << likeAChecksum[i] << " != readChecksum " << readChecksum[i] << " "
464  << " likeAChecksum2 " << likeAChecksum2[i] << " != readChecksum2 " << readChecksum2[i] << ")\n";
465  }
466  }
467 
468  m_residuals_ok.resize(numResiduals(), true);
469 
470  delete[] likeAChecksum;
471  delete[] likeAChecksum2;
472 }
473 
475  double window = 100.;
476  if (which == 0)
477  window = 2. * 30.;
478  else if (which == 1)
479  window = 2. * 30.;
480  else if (which == 2)
481  window = 2. * 20.;
482  else if (which == 3)
483  window = 2. * 50.;
484 
485  TH1F *hist = dir->make<TH1F>(name.c_str(), "", 200, -window, window);
486  for (std::vector<double *>::const_iterator r = residuals_begin(); r != residuals_end(); ++r)
487  hist->Fill(multiplier * (*r)[which]);
488 }
489 
491  std::string name, TFileDirectory *dir, int which, int whichredchi2, double multiplier) {
492  double window = 100.;
493  if (which == 0)
494  window = 2. * 30.;
495  else if (which == 1)
496  window = 2. * 30.;
497  else if (which == 2)
498  window = 2. * 20.;
499  else if (which == 3)
500  window = 2. * 50.;
501 
502  TH1F *hist = dir->make<TH1F>(name.c_str(), "", 200, -window, window);
503  for (std::vector<double *>::const_iterator r = residuals_begin(); r != residuals_end(); ++r) {
504  double weight = 1. / (*r)[whichredchi2];
505  if (TMath::Prob(1. / weight * 12, 12) < 0.99)
506  hist->Fill(multiplier * (*r)[which], weight);
507  }
508 }
509 
511  // first, make a numeric array while discarding some crazy outliers
512  double *data = new double[numResiduals()];
513  int n = 0;
514  for (std::vector<double *>::const_iterator r = m_residuals.begin(); r != m_residuals.end(); r++)
515  if (fabs((*r)[which]) < 50.) {
516  data[n] = (*r)[which];
517  n++;
518  }
519 
520  // compute "3 normal sigma" and regular interquantile ranges
521  const int n_quantiles = 7;
522  double probabilities[n_quantiles] = {0.00135, 0.02275, 0.25, 0.5, 0.75, 0.97725, 0.99865}; // "3 normal sigma"
523  //double probabilities[n_quantiles] = {0.02275, 0.25, 0.75, 0.97725}; // "2 normal sigma"
524  double quantiles[n_quantiles];
525  std::sort(data, data + n);
526  TMath::Quantiles(n, n_quantiles, data, quantiles, probabilities, true, nullptr, 7);
527  delete[] data;
528  double iqr = quantiles[4] - quantiles[2];
529 
530  // estimate optimal bin size according to Freedman-Diaconis rule
531  double hbin = 2 * iqr / pow(n, 1. / 3);
532 
533  a = quantiles[1];
534  b = quantiles[5];
535  nbins = (int)((b - a) / hbin + 3.); // add extra safety margin of 3
536 
537  std::cout << " quantiles: ";
538  for (int i = 0; i < n_quantiles; i++)
539  std::cout << quantiles[i] << " ";
540  std::cout << std::endl;
541  //cout<<"n="<<select_count<<" quantiles ["<<quantiles[1]<<", "<<quantiles[2]<<"] IQR="<<iqr
542  // <<" full range=["<<minx<<","<<maxx<<"]"<<" 2 normal sigma quantile range = ["<<quantiles[0]<<", "<<quantiles[3]<<"]"<<endl;
543  std::cout << " optimal h=" << hbin << " nbins=" << nbins << std::endl;
544 }
545 
546 void MuonResidualsFitter::histogramChi2GaussianFit(int which, double &fit_mean, double &fit_sigma) {
547  int nbins;
548  double a, b;
550  if (a == b || a > b) {
551  fit_mean = a;
552  fit_sigma = 0;
553  return;
554  }
555 
556  TH1D *hist = new TH1D("htmp", "", nbins, a, b);
557  for (std::vector<double *>::const_iterator r = m_residuals.begin(); r != m_residuals.end(); ++r)
558  hist->Fill((*r)[which]);
559 
560  // do simple chi2 gaussian fit
561  TF1 *f1 = new TF1("f1", "gaus", a, b);
562  f1->SetParameter(0, hist->GetEntries());
563  f1->SetParameter(1, 0);
564  f1->SetParameter(2, hist->GetRMS());
565  hist->Fit("f1", "RQ");
566  // hist->Fit(f1,"RQ");
567 
568  fit_mean = f1->GetParameter(1);
569  fit_sigma = f1->GetParameter(2);
570  std::cout << " h(" << nbins << "," << a << "," << b << ") mu=" << fit_mean << " sig=" << fit_sigma << std::endl;
571 
572  delete f1;
573  delete hist;
574 }
575 
576 // simple non-turned ellipsoid selection
577 // THIS WAS selectPeakResiduals_simple, but I changed it to use this simple function as default
578 // void MuonResidualsFitter::selectPeakResiduals_simple(double nsigma, int nvar, int *vars)
579 void MuonResidualsFitter::selectPeakResiduals(double nsigma, int nvar, int *vars) {
580  // does not make sense for small statistics
581  if (numResiduals() < 25)
582  return;
583 
584  int nbefore = numResiduals();
585 
586  //just to be sure (can't see why it might ever be more then 10)
587  assert(nvar <= 10);
588 
589  // estimate nvar-D ellipsoid center & axes
590  for (int v = 0; v < nvar; v++) {
591  int which = vars[v];
593  m_radii[which] = nsigma * m_radii[which];
594  }
595 
596  // filter out residuals that don't fit into the ellipsoid
597  std::vector<double *>::iterator r = m_residuals.begin();
598  while (r != m_residuals.end()) {
599  double ellipsoid_sum = 0;
600  for (int v = 0; v < nvar; v++) {
601  int which = vars[v];
602  if (m_radii[which] == 0.)
603  continue;
604  ellipsoid_sum += pow(((*r)[which] - m_center[which]) / m_radii[which], 2);
605  }
606  if (ellipsoid_sum <= 1.)
607  ++r;
608  else {
609  m_residuals_ok[r - m_residuals.begin()] = false;
610  ++r;
611  // delete [] (*r);
612  // r = m_residuals.erase(r);
613  }
614  }
615  std::cout << " N residuals " << nbefore << " -> " << numResiduals() << std::endl;
616 }
617 
618 // pre-selection using robust covariance estimator
619 // THIS WAS selectPeakResiduals but I changed it to use OTHER simple function as default
620 void MuonResidualsFitter::selectPeakResiduals_simple(double nsigma, int nvar, int *vars) {
621  //std::cout<<"doing selectpeakresiduals: nsig="<<nsigma<<" nvar="<<nvar<<" vars=";
622  for (int i = 0; i < nvar; ++i)
623  std::cout << vars[i] << " ";
624  std::cout << std::endl;
625 
626  // does not make sense for small statistics set to 50
627  // YP changed it to 10 for test
628  if (numResiduals() < 10)
629  return;
630  // if (numResiduals()<50) return;
631 
632  size_t nbefore = numResiduals();
633  std::cout << " N residuals " << nbefore << " ~ "
634  << (size_t)std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true) << std::endl;
635  //just to be sure (can't see why it might ever be more then 10)
636  assert(nvar <= 10 && nvar > 0);
637 
638  std::vector<double *>::iterator r = m_residuals.begin();
639 
640  // it's awkward, but the 1D case has to be handled separately
641  if (nvar == 1) {
642  std::cout << "1D case" << std::endl;
643  // get robust estimates for the peak and sigma
644  double *data = new double[nbefore];
645  for (size_t i = 0; i < nbefore; i++)
646  data[i] = m_residuals[i][vars[0]];
647  double peak, sigma;
648  TRobustEstimator re;
649  re.EvaluateUni(nbefore, data, peak, sigma);
650 
651  // filter out residuals that are more then nsigma away from the peak
652  while (r != m_residuals.end()) {
653  double distance = fabs(((*r)[vars[0]] - peak) / sigma);
654  if (distance <= nsigma)
655  ++r;
656  else {
657  m_residuals_ok[r - m_residuals.begin()] = false;
658  ++r;
659  //delete [] (*r);
660  //r = m_residuals.erase(r);
661  }
662  }
663  std::cout << " N residuals " << nbefore << " -> " << numResiduals() << std::endl;
664  return;
665  } // end 1D case
666 
667  // initialize and run the robust estimator for D>1
668  std::cout << "D>1 case" << std::endl;
669  TRobustEstimator re(nbefore + 1, nvar);
670  std::cout << "nbefore " << nbefore << " nvar " << nvar << std::endl;
671  r = m_residuals.begin();
672  std::cout << "+++++ JUST before loop while (r != m_residuals.end())" << std::endl;
673  int counter1 = 0;
674  while (r != m_residuals.end()) {
675  double *row = new double[nvar];
676  for (int v = 0; v < nvar; v++)
677  row[v] = (*r)[vars[v]];
678  re.AddRow(row);
679  // delete[] row;
680  ++r;
681  counter1++;
682  }
683  std::cout << "counter1 " << counter1 << std::endl;
684  std::cout << "+++++ JUST after loop while (r != m_residuals.end())" << std::endl;
685  re.Evaluate();
686  std::cout << "+++++ JUST after re.Evaluate()" << std::endl;
687 
688  // get nvar-dimensional ellipsoid center & covariance
689  TVectorD M(nvar);
690  re.GetMean(M);
691  TMatrixDSym Cov(nvar);
692  re.GetCovariance(Cov);
693  Cov.Invert();
694 
695  // calculate the normalized radius for this nvar-dimensional ellipsoid from a 1D-Gaussian nsigma equivalent distance
696  double conf_1d = TMath::Erf(nsigma / sqrt(2));
697  double surf_radius = sqrt(TMath::ChisquareQuantile(conf_1d, nvar));
698 
699  // filter out residuals that are outside of the covariance ellipsoid with the above normalized radius
700  r = m_residuals.begin();
701  while (r != m_residuals.end()) {
702  TVectorD res(nvar);
703  for (int v = 0; v < nvar; v++)
704  res[v] = (*r)[vars[v]];
705  double distance = sqrt(Cov.Similarity(res - M));
706  if (distance <= surf_radius)
707  ++r;
708  else {
709  m_residuals_ok[r - m_residuals.begin()] = false;
710  ++r;
711  //delete [] (*r);
712  //r = m_residuals.erase(r);
713  }
714  }
715  std::cout << " N residuals " << nbefore << " -> "
716  << (size_t)std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true) << std::endl;
717 }
718 
719 void MuonResidualsFitter::fiducialCuts(double xMin, double xMax, double yMin, double yMax, bool fidcut1) {
720  int iResidual = -1;
721 
722  int n_station = 9999;
723  int n_wheel = 9999;
724  int n_sector = 9999;
725 
726  double positionX = 9999.;
727  double positionY = 9999.;
728 
729  double chambw = 9999.;
730  double chambl = 9999.;
731 
732  for (std::vector<double *>::const_iterator r = residuals_begin(); r != residuals_end(); ++r) {
733  iResidual++;
734  if (!m_residuals_ok[iResidual])
735  continue;
736 
737  if ((*r)[15] >
738  0.0001) { // this value is greater than zero (chamber width) for 6DOFs stations 1,2,3 better to change for type()!!!
739  n_station = (*r)[12];
740  n_wheel = (*r)[13];
741  n_sector = (*r)[14];
742  positionX = (*r)[4];
743  positionY = (*r)[5];
744  chambw = (*r)[15];
745  chambl = (*r)[16];
746  } else { // in case of 5DOF residual the residual object index is different
747  n_station = (*r)[10];
748  n_wheel = (*r)[11];
749  n_sector = (*r)[12];
750  positionX = (*r)[2];
751  positionY = (*r)[3];
752  chambw = (*r)[13];
753  chambl = (*r)[14];
754  }
755 
756  if (fidcut1) { // this is the standard fiducial cut used so far 80x80 cm in x,y
757  if (positionX >= xMax || positionX <= xMin)
758  m_residuals_ok[iResidual] = false;
759  if (positionY >= yMax || positionY <= yMin)
760  m_residuals_ok[iResidual] = false;
761  }
762 
763  // Implementation of new fiducial cut
764 
765  double dtrkchamx = (chambw / 2.) - positionX; // variables to cut tracks on the edge of the chambers
766  double dtrkchamy = (chambl / 2.) - positionY;
767 
768  if (!fidcut1) {
769  if (n_station == 4) {
770  if ((n_wheel == -1 && n_sector == 3) ||
771  (n_wheel == 1 &&
772  n_sector == 4)) { // FOR SHORT CHAMBER LENGTH IN: WHEEL 1 SECTOR 4 AND WHEEL -1 SECTOR 3
773  if ((n_sector == 1 || n_sector == 2 || n_sector == 3 || n_sector == 5 || n_sector == 6 || n_sector == 7 ||
774  n_sector == 8 || n_sector == 12) &&
775  ((dtrkchamx < 40 || dtrkchamx > 380) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
776  m_residuals_ok[iResidual] = false;
777  if ((n_sector == 4 || n_sector == 13) &&
778  ((dtrkchamx < 40 || dtrkchamx > 280) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
779  m_residuals_ok[iResidual] = false;
780  if ((n_sector == 9 || n_sector == 11) &&
781  ((dtrkchamx < 40 || dtrkchamx > 180) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
782  m_residuals_ok[iResidual] = false;
783  if ((n_sector == 10 || n_sector == 14) &&
784  ((dtrkchamx < 40 || dtrkchamx > 220) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
785  m_residuals_ok[iResidual] = false;
786  } else {
787  if ((n_sector == 1 || n_sector == 2 || n_sector == 3 || n_sector == 5 || n_sector == 6 || n_sector == 7 ||
788  n_sector == 8 || n_sector == 12) &&
789  ((dtrkchamx < 40 || dtrkchamx > 380) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
790  m_residuals_ok[iResidual] = false;
791  if ((n_sector == 4 || n_sector == 13) &&
792  ((dtrkchamx < 40 || dtrkchamx > 280) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
793  m_residuals_ok[iResidual] = false;
794  if ((n_sector == 9 || n_sector == 11) &&
795  ((dtrkchamx < 40 || dtrkchamx > 180) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
796  m_residuals_ok[iResidual] = false;
797  if ((n_sector == 10 || n_sector == 14) &&
798  ((dtrkchamx < 40 || dtrkchamx > 220) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
799  m_residuals_ok[iResidual] = false;
800  }
801  } else {
802  if ((n_wheel == -1 && n_sector == 3) || (n_wheel == 1 && n_sector == 4)) {
803  if (n_station == 1 && ((dtrkchamx < 30.0 || dtrkchamx > 190.0) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
804  m_residuals_ok[iResidual] = false;
805  if (n_station == 2 && ((dtrkchamx < 30.0 || dtrkchamx > 240.0) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
806  m_residuals_ok[iResidual] = false;
807  if (n_station == 3 && ((dtrkchamx < 30.0 || dtrkchamx > 280.0) || (dtrkchamy < 40.0 || dtrkchamy > 170.0)))
808  m_residuals_ok[iResidual] = false;
809  } else {
810  if (n_station == 1 && ((dtrkchamx < 30.0 || dtrkchamx > 190.0) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
811  m_residuals_ok[iResidual] = false;
812  if (n_station == 2 && ((dtrkchamx < 30.0 || dtrkchamx > 240.0) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
813  m_residuals_ok[iResidual] = false;
814  if (n_station == 3 && ((dtrkchamx < 30.0 || dtrkchamx > 280.0) || (dtrkchamy < 40.0 || dtrkchamy > 210.0)))
815  m_residuals_ok[iResidual] = false;
816  }
817  }
818  }
819  }
820 }
821 
822 void MuonResidualsFitter::correctBField(int idx_momentum, int idx_q) {
823  const int Nbin = 17;
824  // find max 1/pt and bin width
825  double min_pt = 9999.;
826  for (std::vector<double *>::const_iterator r = residuals_begin(); r != residuals_end(); ++r) {
827  double pt = fabs((*r)[idx_momentum]);
828  if (pt < min_pt)
829  min_pt = pt;
830  }
831  min_pt -= 0.01; // to prevent bin # overflow
832  const double bin_width = 1. / min_pt / Nbin;
833 
834  // fill indices of positive and negative charge residuals in each bin
835  std::vector<size_t> pos[Nbin], neg[Nbin], to_erase;
836  for (size_t i = 0; i < m_residuals.size(); i++) {
837  if (!m_residuals_ok[i])
838  continue;
839  int bin = (int)floor(1. / fabs(m_residuals[i][idx_momentum]) / bin_width);
840  if (m_residuals[i][idx_q] > 0)
841  pos[bin].push_back(i);
842  else
843  neg[bin].push_back(i);
844  }
845 
846  // equalize pos and neg in each bin
847  for (int j = 0; j < Nbin; j++) {
848  size_t psize = pos[j].size();
849  size_t nsize = neg[j].size();
850  if (psize == nsize)
851  continue;
852 
853  std::set<int> idx_set; // use a set to collect certain number of unique random indices to erase
854  if (psize > nsize) {
855  while (idx_set.size() < psize - nsize)
856  idx_set.insert(gRandom->Integer(psize));
857  for (std::set<int>::iterator it = idx_set.begin(); it != idx_set.end(); it++)
858  to_erase.push_back(pos[j][*it]);
859  } else {
860  while (idx_set.size() < nsize - psize)
861  idx_set.insert(gRandom->Integer(nsize));
862  for (std::set<int>::iterator it = idx_set.begin(); it != idx_set.end(); it++)
863  to_erase.push_back(neg[j][*it]);
864  }
865  }
866  // sort in descending order, so we safely go from higher to lower indices:
867  std::sort(to_erase.begin(), to_erase.end(), std::greater<int>());
868  for (std::vector<size_t>::const_iterator e = to_erase.begin(); e != to_erase.end(); ++e) {
869  m_residuals_ok[*e] = false;
870  //delete[] *(m_residuals.begin() + *e);
871  //m_residuals.erase(m_residuals.begin() + *e);
872  }
873 
874  std::vector<size_t> apos[Nbin], aneg[Nbin];
875  for (size_t i = 0; i < m_residuals.size(); i++) {
876  if (!m_residuals_ok[i])
877  continue;
878  int bin = (int)floor(1. / fabs(m_residuals[i][idx_momentum]) / bin_width);
879  if (m_residuals[i][idx_q] > 0)
880  apos[bin].push_back(i);
881  else
882  aneg[bin].push_back(i);
883  }
884  for (int j = 0; j < Nbin; j++)
885  std::cout << "bin " << j << ": [pos,neg] sizes before & after: [" << pos[j].size() << "," << neg[j].size()
886  << "] -> [" << apos[j].size() << "," << aneg[j].size() << "]" << std::endl;
887  std::cout << " N residuals " << m_residuals.size() << " -> "
888  << (size_t)std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true) << std::endl;
889 }
890 
892  // it should probably be faster then doing erase
893  size_t n_ok = (size_t)std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true);
894  std::vector<double *> tmp(n_ok, nullptr);
895  std::cout << "residuals sizes: all=" << m_residuals.size() << " good=" << n_ok << std::endl;
896  int iok = 0;
897  for (size_t i = 0; i < m_residuals.size(); i++) {
898  if (!m_residuals_ok[i])
899  continue;
900  tmp[iok++] = m_residuals[i];
901  }
902  m_residuals.swap(tmp);
903 
904  std::vector<bool> tmp_ok(n_ok, true);
905  m_residuals_ok.swap(tmp_ok);
906 
907  std::cout << "residuals size after eraseNotSelectedResiduals =" << m_residuals.size()
908  << " ok size=" << m_residuals_ok.size() << std::endl;
909 }
MuonResidualsFitter::fiducialCuts
void fiducialCuts(double xMin=-80.0, double xMax=80.0, double yMin=-80.0, double yMax=80.0, bool fidcut1=false)
Definition: MuonResidualsFitter.cc:719
MuonResidualsFitter::correctBField
virtual void correctBField()=0
MuonResidualsFitter_TMinuit
static TMinuit * MuonResidualsFitter_TMinuit
Definition: MuonResidualsFitter.cc:28
mps_fire.i
i
Definition: mps_fire.py:428
start
Definition: start.py:1
MuonResidualsFitter_pureGaussian_TF1
Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:38
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
MuonResidualsFitter::MuonResidualsFitter
MuonResidualsFitter(int residualsModel, int minHits, int useResiduals, bool weightAlignment=true)
Definition: MuonResidualsFitter.cc:151
MuonResidualsFitter::npar
virtual int npar()=0
align_cfg.useResiduals
useResiduals
Definition: align_cfg.py:36
MuonResidualsFitter::m_residuals_ok
std::vector< bool > m_residuals_ok
Definition: MuonResidualsFitter.h:233
step
step
Definition: StallMonitor.cc:94
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
MuonResidualsFitter::dofit
bool dofit(void(*fcn)(int &, double *, double &, double *, int), std::vector< int > &parNum, std::vector< std::string > &parName, std::vector< double > &start, std::vector< double > &step, std::vector< double > &low, std::vector< double > &high)
Definition: MuonResidualsFitter.cc:272
gather_cfg.cout
cout
Definition: gather_cfg.py:144
pos
Definition: PixelAliasList.h:18
fitWZ.arglist
arglist
Definition: fitWZ.py:39
MuonResidualsFitter::m_fixed
std::vector< bool > m_fixed
Definition: MuonResidualsFitter.h:228
MuonResidualsFitter::~MuonResidualsFitter
virtual ~MuonResidualsFitter()
Definition: MuonResidualsFitter.cc:165
svgfig.window
def window(xmin, xmax, ymin, ymax, x=0, y=0, width=100, height=100, xlogbase=None, ylogbase=None, minusInfinity=-1000, flipx=False, flipy=True)
Definition: svgfig.py:643
MuonResidualsFitter::numResiduals
long numResiduals() const
Definition: MuonResidualsFitter.h:125
cms::cuda::assert
assert(be >=bs)
CustomPhysics_cfi.gamma
gamma
Definition: CustomPhysics_cfi.py:17
MuonResidualsFitter::eraseNotSelectedResiduals
void eraseNotSelectedResiduals()
Definition: MuonResidualsFitter.cc:891
MuonResidualsFitter::kGaussPowerTails
Definition: MuonResidualsFitter.h:86
MuonResidualsFitter::parNum2parIdx
int parNum2parIdx(int parNum)
Definition: MuonResidualsFitter.h:156
photonAnalyzer_cfi.yMin
yMin
Definition: photonAnalyzer_cfi.py:86
TFileDirectory
Definition: TFileDirectory.h:24
fcn
void fcn(int &, double *, double &, double *, int)
Definition: LASBarrelAlgorithm.cc:388
findQualityFiles.v
v
Definition: findQualityFiles.py:179
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
tools.TF1
TF1
Definition: tools.py:23
cosmictrackSelector_cfi.min_pt
min_pt
Definition: cosmictrackSelector_cfi.py:18
MuonResidualsFitter_logGaussPowerTails
double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma)
Definition: MuonResidualsFitter.cc:128
align_cfg.residualsModel
residualsModel
Definition: align_cfg.py:34
photonAnalyzer_cfi.xMin
xMin
Definition: photonAnalyzer_cfi.py:82
MuonResidualsFitter_compute_log_convolution
double MuonResidualsFitter_compute_log_convolution(double toversigma, double gammaoversigma, double max, double step, double power)
Definition: MuonResidualsFitter.cc:54
MuonResidualsFitter::plotweighted
void plotweighted(std::string name, TFileDirectory *dir, int which, int whichredchi2, double multiplier)
Definition: MuonResidualsFitter.cc:490
MuonResidualsFitter_lookup_table
double MuonResidualsFitter_lookup_table[MuonResidualsFitter_numgsbins][MuonResidualsFitter_numtsbins]
Definition: MuonResidualsFitter.cc:26
alignCSCRings.s
s
Definition: alignCSCRings.py:92
reco::ceil
constexpr int32_t ceil(float num)
Definition: constexpr_cmath.h:7
compare.hist
hist
Definition: compare.py:376
vars
vars
Definition: DeepTauId.cc:159
MuonResidualsFitterFitInfo
Definition: MuonResidualsFitter.h:250
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
MuonResidualsFitter::kPowerLawTails
Definition: MuonResidualsFitter.h:86
MuonResidualsFitter::selectPeakResiduals_simple
void selectPeakResiduals_simple(double nsigma, int nvar, int *vars)
Definition: MuonResidualsFitter.cc:620
LaserClient_cfi.nbins
nbins
Definition: LaserClient_cfi.py:51
submitPVResolutionJobs.count
count
Definition: submitPVResolutionJobs.py:352
MuonResidualsFitter::m_residualsModel
int m_residualsModel
Definition: MuonResidualsFitter.h:224
MuonResidualsFitter.h
MuonResidualsFitter_ROOTVoigt_TF1
Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:124
multiplicitycorr_cfi.yMax
yMax
Definition: multiplicitycorr_cfi.py:6
b
double b
Definition: hdecay.h:118
MuonResidualsFitter::ndata
virtual int ndata()=0
MuonResidualsFitter::write
void write(FILE *file, int which=0)
Definition: MuonResidualsFitter.cc:383
MuonResidualsFitter_numgsbins
const int MuonResidualsFitter_numgsbins
Definition: MuonResidualsFitter.cc:22
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
MuonResidualsFitter::covarianceElement
double covarianceElement(int parNum1, int parNum2)
Definition: MuonResidualsFitter.cc:191
MuonResidualsFitter_table_initialized
bool MuonResidualsFitter_table_initialized
Definition: MuonResidualsFitter.cc:25
MuonResidualsFitter::kPureGaussian
Definition: MuonResidualsFitter.h:86
a
double a
Definition: hdecay.h:119
MuonResidualsFitter::plotsimple
void plotsimple(std::string name, TFileDirectory *dir, int which, double multiplier)
Definition: MuonResidualsFitter.cc:474
MuonResidualsFitter_gsbinsize
const double MuonResidualsFitter_gsbinsize
Definition: MuonResidualsFitter.cc:20
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
MuonResidualsFitter::m_residuals
std::vector< double * > m_residuals
Definition: MuonResidualsFitter.h:232
MuonResidualsFitter::residuals_end
std::vector< double * >::const_iterator residuals_end() const
Definition: MuonResidualsFitter.h:191
MuonResidualsFitter::m_error
std::vector< double > m_error
Definition: MuonResidualsFitter.h:236
createfilelist.int
int
Definition: createfilelist.py:10
FrontierConditions_GlobalTag_cff.file
file
Definition: FrontierConditions_GlobalTag_cff.py:13
MuonResidualsFitter_powerLawTails_TF1
Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:116
fftjetcommon_cfi.sy
sy
Definition: fftjetcommon_cfi.py:203
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
MuonResidualsFitter_GaussPowerTails_TF1
Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par)
Definition: MuonResidualsFitter.cc:142
MuonResidualsFitter_logPureGaussian
double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma)
Definition: MuonResidualsFitter.cc:31
PVValHelper::dy
Definition: PVValidationHelpers.h:49
MuonResidualsFitter::computeHistogramRangeAndBinning
void computeHistogramRangeAndBinning(int which, int &nbins, double &a, double &b)
Definition: MuonResidualsFitter.cc:510
MuonResidualsFitter::m_strategy
int m_strategy
Definition: MuonResidualsFitter.h:230
MuonResidualsFitter::fill
void fill(double *residual)
Definition: MuonResidualsFitter.cc:186
res
Definition: Electron.h:6
MuonResidualsFitter::histogramChi2GaussianFit
void histogramChi2GaussianFit(int which, double &fit_mean, double &fit_sigma)
Definition: MuonResidualsFitter.cc:546
MuonResidualsFitter::kROOTVoigt
Definition: MuonResidualsFitter.h:86
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:34
MuonResidualsFitter_numtsbins
const int MuonResidualsFitter_numtsbins
Definition: MuonResidualsFitter.cc:23
alignCSCRings.r
r
Definition: alignCSCRings.py:93
newFWLiteAna.bin
bin
Definition: newFWLiteAna.py:161
MuonResidualsFitter::m_loglikelihood
double m_loglikelihood
Definition: MuonResidualsFitter.h:238
MuonResidualsFitter::m_value
std::vector< double > m_value
Definition: MuonResidualsFitter.h:235
LaserClient_cfi.high
high
Definition: LaserClient_cfi.py:50
heppy_batch.val
val
Definition: heppy_batch.py:351
MuonResidualsFitter_logROOTVoigt
double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma)
Definition: MuonResidualsFitter.cc:120
Exception
Definition: hltDiff.cc:246
align_cfg.weightAlignment
tuple weightAlignment
Definition: align_cfg.py:30
genVertex_cff.x
x
Definition: genVertex_cff.py:12
MuonResidualsFitter_integrate_pureGaussian
double MuonResidualsFitter_integrate_pureGaussian(double low, double high, double center, double sigma)
Definition: MuonResidualsFitter.cc:146
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
MuonResidualsFitter::fix
void fix(int parNum, bool dofix=true)
Definition: MuonResidualsFitter.cc:171
MuonResidualsFitter::selectPeakResiduals
void selectPeakResiduals(double nsigma, int nvar, int *vars)
Definition: MuonResidualsFitter.cc:579
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
multiplicitycorr_cfi.xMax
xMax
Definition: multiplicitycorr_cfi.py:5
MuonResidualsFitter::read
void read(FILE *file, int which=0)
Definition: MuonResidualsFitter.cc:418
MuonResidualsFitter::initialize_table
void initialize_table()
Definition: MuonResidualsFitter.cc:205
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
data
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
MuonResidualsFitter::m_center
double m_center[20]
Definition: MuonResidualsFitter.h:245
MuonResidualsFitter_tsbinsize
const double MuonResidualsFitter_tsbinsize
Definition: MuonResidualsFitter.cc:21
MuonResidualsFitter_logPowerLawTails
double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma)
Definition: MuonResidualsFitter.cc:85
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
MuonResidualsFitter_logPureGaussian2D
double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r)
Definition: MuonResidualsFitter.cc:43
cms::Exception
Definition: Exception.h:70
MuonResidualsFitter::kPureGaussian2D
Definition: MuonResidualsFitter.h:86
MuonResidualsFitter::m_cov
TMatrixDSym m_cov
Definition: MuonResidualsFitter.h:237
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
MuonResidualsFitter::residualsModel
int residualsModel() const
Definition: MuonResidualsFitter.h:124
MuonResidualsFitter::m_radii
double m_radii[20]
Definition: MuonResidualsFitter.h:246
MuonResidualsFitter::fixed
bool fixed(int parNum)
Definition: MuonResidualsFitter.cc:178
DeadROC_duringRun.f1
f1
Definition: DeadROC_duringRun.py:219
postprocess-scan-build.rows
rows
Definition: postprocess-scan-build.py:11
HLT_FULL_cff.distance
distance
Definition: HLT_FULL_cff.py:7799
MuonResidualsFitter::numsegments
long numsegments()
Definition: MuonResidualsFitter.cc:198
eostools.which
def which(cmd)
Definition: eostools.py:336
MuonResidualsFitter::residuals_begin
std::vector< double * >::const_iterator residuals_begin() const
Definition: MuonResidualsFitter.h:190
LaserClient_cfi.low
low
Definition: LaserClient_cfi.py:52
PVValHelper::dx
Definition: PVValidationHelpers.h:48
weight
Definition: weight.py:1
reco_skim_cfg_mod.minHits
minHits
Definition: reco_skim_cfg_mod.py:177
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
MuonResidualsFitter::inform
virtual void inform(TMinuit *tMinuit)=0
fftjetcommon_cfi.sx
sx
Definition: fftjetcommon_cfi.py:202