CMS 3D CMS Logo

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