00001
00002
00003 #ifdef STANDALONE_FITTER
00004 #include "MuonResidualsFitter.h"
00005 #else
00006 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsFitter.h"
00007 #endif
00008
00009 #include <fstream>
00010 #include <set>
00011 #include "TMath.h"
00012 #include "TH1.h"
00013 #include "TF1.h"
00014 #include "TRobustEstimator.h"
00015
00016
00017
00018
00019
00020 const double MuonResidualsFitter_gsbinsize = 0.01;
00021 const double MuonResidualsFitter_tsbinsize = 0.1;
00022 const int MuonResidualsFitter_numgsbins = 500;
00023 const int MuonResidualsFitter_numtsbins = 500;
00024
00025 bool MuonResidualsFitter_table_initialized = false;
00026 double MuonResidualsFitter_lookup_table[MuonResidualsFitter_numgsbins][MuonResidualsFitter_numtsbins];
00027
00028 static TMinuit* MuonResidualsFitter_TMinuit;
00029
00030
00031 double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma)
00032 {
00033 sigma = fabs(sigma);
00034 static const double cgaus = 0.5 * log( 2.*M_PI );
00035 return (-pow(residual - center, 2) *0.5 / sigma / sigma) - cgaus - log(sigma);
00036 }
00037
00038
00039 Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par)
00040 {
00041 return par[0] * exp(MuonResidualsFitter_logPureGaussian(xvec[0], par[1], par[2]));
00042 }
00043
00044
00045
00046 double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r)
00047 {
00048 sx = fabs(sx); sy = fabs(sy);
00049 static const double c2gaus = log( 2.*M_PI );
00050 double one_r2 = 1. - r*r;
00051 double dx = x - x0;
00052 double dy = y - y0;
00053 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));
00054 }
00055
00056
00057 double MuonResidualsFitter_compute_log_convolution(double toversigma, double gammaoversigma, double max, double step, double power)
00058 {
00059 if (gammaoversigma == 0.) return (-toversigma*toversigma/2.) - log(sqrt(2*M_PI));
00060
00061 double sum = 0.;
00062 double uplus = 0.;
00063 double integrandplus_last = 0.;
00064 double integrandminus_last = 0.;
00065 for (double inc = 0.; uplus < max; inc += step)
00066 {
00067 double uplus_last = uplus;
00068 uplus = pow(inc, power);
00069
00070 double integrandplus = exp(-pow(uplus - toversigma, 2)/2.) / (uplus*uplus/gammaoversigma + gammaoversigma) / M_PI / sqrt(2.*M_PI);
00071 double integrandminus = exp(-pow(-uplus - toversigma, 2)/2.) / (uplus*uplus/gammaoversigma + gammaoversigma) / M_PI / sqrt(2.*M_PI);
00072
00073 sum += integrandplus * (uplus - uplus_last);
00074 sum += 0.5 * fabs(integrandplus_last - integrandplus) * (uplus - uplus_last);
00075
00076 sum += integrandminus * (uplus - uplus_last);
00077 sum += 0.5 * fabs(integrandminus_last - integrandminus) * (uplus - uplus_last);
00078
00079 integrandplus_last = integrandplus;
00080 integrandminus_last = integrandminus;
00081 }
00082 return log(sum);
00083 }
00084
00085
00086 double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma)
00087 {
00088 sigma = fabs(sigma);
00089 double gammaoversigma = fabs(gamma / sigma);
00090 double toversigma = fabs((residual - center) / sigma);
00091
00092 int gsbin0 = int(floor(gammaoversigma / MuonResidualsFitter_gsbinsize));
00093 int gsbin1 = int(ceil(gammaoversigma / MuonResidualsFitter_gsbinsize));
00094 int tsbin0 = int(floor(toversigma / MuonResidualsFitter_tsbinsize));
00095 int tsbin1 = int(ceil(toversigma / MuonResidualsFitter_tsbinsize));
00096
00097 bool gsisbad = (gsbin0 >= MuonResidualsFitter_numgsbins || gsbin1 >= MuonResidualsFitter_numgsbins);
00098 bool tsisbad = (tsbin0 >= MuonResidualsFitter_numtsbins || tsbin1 >= MuonResidualsFitter_numtsbins);
00099
00100 if (gsisbad || tsisbad)
00101 {
00102 return log(gammaoversigma/M_PI) - log(toversigma*toversigma + gammaoversigma*gammaoversigma) - log(sigma);
00103 }
00104 else
00105 {
00106 double val00 = MuonResidualsFitter_lookup_table[gsbin0][tsbin0];
00107 double val01 = MuonResidualsFitter_lookup_table[gsbin0][tsbin1];
00108 double val10 = MuonResidualsFitter_lookup_table[gsbin1][tsbin0];
00109 double val11 = MuonResidualsFitter_lookup_table[gsbin1][tsbin1];
00110
00111 double val0 = val00 + ((toversigma / MuonResidualsFitter_tsbinsize) - tsbin0) * (val01 - val00);
00112 double val1 = val10 + ((toversigma / MuonResidualsFitter_tsbinsize) - tsbin0) * (val11 - val10);
00113
00114 double val = val0 + ((gammaoversigma / MuonResidualsFitter_gsbinsize) - gsbin0) * (val1 - val0);
00115
00116 return val - log(sigma);
00117 }
00118 }
00119
00120
00121
00122 Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par)
00123 {
00124 return par[0] * exp(MuonResidualsFitter_logPowerLawTails(xvec[0], par[1], par[2], par[3]));
00125 }
00126
00127
00128 double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma)
00129 {
00130 return log(TMath::Voigt(residual - center, fabs(sigma), fabs(gamma)*2.));
00131 }
00132
00133
00134 Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par)
00135 {
00136 return par[0] * TMath::Voigt(xvec[0] - par[1], fabs(par[2]), fabs(par[3])*2.);
00137 }
00138
00139
00140 double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma) {
00141 double x = residual-center;
00142 double s = fabs(sigma);
00143 double m = 2.*s;
00144 static const double e2 = exp(-2.), sqerf = sqrt(2.*M_PI)*erf(sqrt(2.));
00145 double a = pow(m,4)*e2;
00146 double n = sqerf*s + 2.*a*pow(m,-3)/3.;
00147
00148 if (fabs(x)<m) return -x*x/(2.*s*s) - log(n);
00149 else return log(a) -4.*log(fabs(x)) - log(n);
00150 }
00151
00152
00153 Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par)
00154 {
00155 return par[0] * exp(MuonResidualsFitter_logGaussPowerTails(xvec[0], par[1], par[2]));
00156 }
00157
00158
00159 double MuonResidualsFitter_integrate_pureGaussian(double low, double high, double center, double sigma)
00160 {
00161 static const double isqr2 = 1./sqrt(2.);
00162 return (erf((high + center) * isqr2 / sigma) - erf((low + center) * isqr2 / sigma)) * exp(0.5/sigma/sigma) * 0.5;
00163 }
00164
00165
00166 void MuonResidualsFitter::initialize_table()
00167 {
00168 if (MuonResidualsFitter_table_initialized || residualsModel() != kPowerLawTails) return;
00169 MuonResidualsFitter_table_initialized = true;
00170
00171 std::ifstream convolution_table("convolution_table.txt");
00172 if (convolution_table.is_open())
00173 {
00174 int numgsbins = 0;
00175 int numtsbins = 0;
00176 double tsbinsize = 0.;
00177 double gsbinsize = 0.;
00178
00179 convolution_table >> numgsbins >> numtsbins >> tsbinsize >> gsbinsize;
00180 if (numgsbins != MuonResidualsFitter_numgsbins || numtsbins != MuonResidualsFitter_numtsbins ||
00181 tsbinsize != MuonResidualsFitter_tsbinsize || gsbinsize != MuonResidualsFitter_gsbinsize)
00182 {
00183 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";
00184 }
00185
00186 for (int gsbin = 0; gsbin < MuonResidualsFitter_numgsbins; gsbin++)
00187 {
00188 for (int tsbin = 0; tsbin < MuonResidualsFitter_numtsbins; tsbin++)
00189 {
00190 int read_gsbin = 0;
00191 int read_tsbin = 0;
00192 double val = 0.;
00193
00194 convolution_table >> read_gsbin >> read_tsbin >> val;
00195 if (read_gsbin != gsbin || read_tsbin != tsbin)
00196 {
00197 throw cms::Exception("MuonResidualsFitter") << "convolution_table.txt is out of order. Throw it away and let the fitter re-create the file.\n";
00198 }
00199 MuonResidualsFitter_lookup_table[gsbin][tsbin] = val;
00200 }
00201 }
00202 convolution_table.close();
00203 }
00204 else
00205 {
00206 std::ofstream convolution_table2("convolution_table.txt");
00207
00208 if (!convolution_table2.is_open()) throw cms::Exception("MuonResidualsFitter") << "Couldn't write to file convolution_table.txt\n";
00209
00210 convolution_table2 << MuonResidualsFitter_numgsbins << " " << MuonResidualsFitter_numtsbins << " " << MuonResidualsFitter_tsbinsize << " " << MuonResidualsFitter_gsbinsize << std::endl;
00211
00212 std::cout << "Initializing convolution look-up table (takes a few minutes)..." << std::endl;
00213
00214 for (int gsbin = 0; gsbin < MuonResidualsFitter_numgsbins; gsbin++)
00215 {
00216 double gammaoversigma = double(gsbin) * MuonResidualsFitter_gsbinsize;
00217 std::cout << " gsbin " << gsbin << "/" << MuonResidualsFitter_numgsbins << std::endl;
00218 for (int tsbin = 0; tsbin < MuonResidualsFitter_numtsbins; tsbin++)
00219 {
00220 double toversigma = double(tsbin) * MuonResidualsFitter_tsbinsize;
00221
00222
00223 MuonResidualsFitter_lookup_table[gsbin][tsbin] = MuonResidualsFitter_compute_log_convolution(toversigma, gammaoversigma);
00224
00225
00226
00227
00228 convolution_table2 << gsbin << " " << tsbin << " " << MuonResidualsFitter_lookup_table[gsbin][tsbin] << std::endl;
00229 }
00230 }
00231 convolution_table2.close();
00232 std::cout << "Initialization done!" << std::endl;
00233 }
00234 }
00235
00236
00237 bool MuonResidualsFitter::dofit(void (*fcn)(int&,double*,double&,double*,int), std::vector<int> &parNum, std::vector<std::string> &parName,
00238 std::vector<double> &start, std::vector<double> &step, std::vector<double> &low, std::vector<double> &high)
00239 {
00240 MuonResidualsFitterFitInfo *fitinfo = new MuonResidualsFitterFitInfo(this);
00241
00242 MuonResidualsFitter_TMinuit = new TMinuit(npar());
00243 MuonResidualsFitter_TMinuit->SetPrintLevel(m_printLevel);
00244 MuonResidualsFitter_TMinuit->SetObjectFit(fitinfo);
00245 MuonResidualsFitter_TMinuit->SetFCN(fcn);
00246 inform(MuonResidualsFitter_TMinuit);
00247
00248 std::vector<int>::const_iterator iNum = parNum.begin();
00249 std::vector<std::string>::const_iterator iName = parName.begin();
00250 std::vector<double>::const_iterator istart = start.begin();
00251 std::vector<double>::const_iterator istep = step.begin();
00252 std::vector<double>::const_iterator ilow = low.begin();
00253 std::vector<double>::const_iterator ihigh = high.begin();
00254
00255
00256
00257 for (; iNum != parNum.end(); ++iNum, ++iName, ++istart, ++istep, ++ilow, ++ihigh)
00258 {
00259 MuonResidualsFitter_TMinuit->DefineParameter(*iNum, iName->c_str(), *istart, *istep, *ilow, *ihigh);
00260 if (fixed(*iNum)) MuonResidualsFitter_TMinuit->FixParameter(*iNum);
00261 }
00262
00263 double arglist[10];
00264 int ierflg;
00265 int smierflg;
00266
00267
00268 for (int i = 0; i < 10; i++) arglist[i] = 0.;
00269 arglist[0] = 0.5;
00270 ierflg = 0;
00271 smierflg = 0;
00272 MuonResidualsFitter_TMinuit->mnexcm("SET ERR", arglist, 1, ierflg);
00273 if (ierflg != 0) { delete MuonResidualsFitter_TMinuit; delete fitinfo; return false; }
00274
00275
00276 for (int i = 0; i < 10; i++) arglist[i] = 0.;
00277 arglist[0] = m_strategy;
00278 ierflg = 0;
00279 MuonResidualsFitter_TMinuit->mnexcm("SET STR", arglist, 1, ierflg);
00280 if (ierflg != 0) { delete MuonResidualsFitter_TMinuit; delete fitinfo; return false; }
00281
00282 bool try_again = false;
00283
00284
00285 for (int i = 0; i < 10; i++) arglist[i] = 0.;
00286 arglist[0] = 50000;
00287 ierflg = 0;
00288 MuonResidualsFitter_TMinuit->mnexcm("MIGRAD", arglist, 1, ierflg);
00289 if (ierflg != 0) try_again = true;
00290
00291
00292 if (try_again)
00293 {
00294 for (int i = 0; i < 10; i++) arglist[i] = 0.;
00295 arglist[0] = 50000;
00296 MuonResidualsFitter_TMinuit->mnexcm("MIGRAD", arglist, 1, smierflg);
00297 }
00298
00299 Double_t fmin, fedm, errdef;
00300 Int_t npari, nparx, istat;
00301 MuonResidualsFitter_TMinuit->mnstat(fmin, fedm, errdef, npari, nparx, istat);
00302
00303 if (istat != 3)
00304 {
00305 for (int i = 0; i < 10; i++) arglist[i] = 0.;
00306 ierflg = 0;
00307 MuonResidualsFitter_TMinuit->mnexcm("HESSE", arglist, 0, ierflg);
00308 }
00309
00310
00311 m_loglikelihood = -fmin;
00312
00313 m_value.clear();
00314 m_error.clear();
00315 for (int i = 0; i < npar(); i++)
00316 {
00317 double v, e;
00318 MuonResidualsFitter_TMinuit->GetParameter(i, v, e);
00319 m_value.push_back(v);
00320 m_error.push_back(e);
00321 }
00322 m_cov.ResizeTo(npar(),npar());
00323 MuonResidualsFitter_TMinuit->mnemat( m_cov.GetMatrixArray(), npar());
00324
00325 delete MuonResidualsFitter_TMinuit;
00326 delete fitinfo;
00327 if (smierflg != 0) return false;
00328 return true;
00329 }
00330
00331
00332 void MuonResidualsFitter::write(FILE *file, int which)
00333 {
00334 long rows = numResiduals();
00335 int cols = ndata();
00336 int whichcopy = which;
00337
00338 fwrite(&rows, sizeof(long), 1, file);
00339 fwrite(&cols, sizeof(int), 1, file);
00340 fwrite(&whichcopy, sizeof(int), 1, file);
00341
00342 double *likeAChecksum = new double[cols];
00343 double *likeAChecksum2 = new double[cols];
00344 for (int i = 0; i < cols; i++)
00345 {
00346 likeAChecksum[i] = 0.;
00347 likeAChecksum2[i] = 0.;
00348 }
00349
00350 for (std::vector<double*>::const_iterator residual = residuals_begin(); residual != residuals_end(); ++residual)
00351 {
00352 fwrite((*residual), sizeof(double), cols, file);
00353 for (int i = 0; i < cols; i++)
00354 {
00355 if (fabs((*residual)[i]) > likeAChecksum[i]) likeAChecksum[i] = fabs((*residual)[i]);
00356 if (fabs((*residual)[i]) < likeAChecksum2[i]) likeAChecksum2[i] = fabs((*residual)[i]);
00357 }
00358 }
00359
00360
00361
00362 fwrite(likeAChecksum, sizeof(double), cols, file);
00363 fwrite(likeAChecksum2, sizeof(double), cols, file);
00364
00365 delete [] likeAChecksum;
00366 delete [] likeAChecksum2;
00367 }
00368
00369 void MuonResidualsFitter::read(FILE *file, int which)
00370 {
00371 long rows = -100;
00372 int cols = -100;
00373 int readwhich = -100;
00374
00375 fread(&rows, sizeof(long), 1, file);
00376 fread(&cols, sizeof(int), 1, file);
00377 fread(&readwhich, sizeof(int), 1, file);
00378
00379 if (cols != ndata() || rows < 0 || readwhich != which)
00380 {
00381 throw cms::Exception("MuonResidualsFitter") << "temporary file is corrupted (which = " << which << " readwhich = " << readwhich << " rows = " << rows << " cols = " << cols << ")\n";
00382 }
00383
00384 double *likeAChecksum = new double[cols];
00385 double *likeAChecksum2 = new double[cols];
00386 for (int i = 0; i < cols; i++)
00387 {
00388 likeAChecksum[i] = 0.;
00389 likeAChecksum2[i] = 0.;
00390 }
00391
00392 m_residuals.reserve(rows);
00393 for (long row = 0; row < rows; row++)
00394 {
00395 double *residual = new double[cols];
00396 fread(residual, sizeof(double), cols, file);
00397 fill(residual);
00398 for (int i = 0; i < cols; i++)
00399 {
00400 if (fabs(residual[i]) > likeAChecksum[i]) likeAChecksum[i] = fabs(residual[i]);
00401 if (fabs(residual[i]) < likeAChecksum2[i]) likeAChecksum2[i] = fabs(residual[i]);
00402 }
00403 }
00404
00405 double *readChecksum = new double[cols];
00406 double *readChecksum2 = new double[cols];
00407 fread(readChecksum, sizeof(double), cols, file);
00408 fread(readChecksum2, sizeof(double), cols, file);
00409
00410 for (int i = 0; i < cols; i++)
00411 {
00412 if (fabs(likeAChecksum[i] - readChecksum[i]) > 1e-10 || fabs(1./likeAChecksum2[i] - 1./readChecksum2[i]) > 1e10)
00413 {
00414 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";
00415 }
00416 }
00417
00418 m_residuals_ok.resize(numResiduals(), true);
00419
00420 delete [] likeAChecksum;
00421 delete [] likeAChecksum2;
00422 }
00423
00424
00425 void MuonResidualsFitter::plotsimple(std::string name, TFileDirectory *dir, int which, double multiplier)
00426 {
00427 double window = 100.;
00428 if (which == 0) window = 2.*30.;
00429 else if (which == 1) window = 2.*30.;
00430 else if (which == 2) window = 2.*20.;
00431 else if (which == 3) window = 2.*50.;
00432
00433 TH1F *hist = dir->make<TH1F>(name.c_str(), "", 200, -window, window);
00434 for (std::vector<double*>::const_iterator r = residuals_begin(); r != residuals_end(); ++r) hist->Fill(multiplier * (*r)[which]);
00435 }
00436
00437
00438 void MuonResidualsFitter::plotweighted(std::string name, TFileDirectory *dir, int which, int whichredchi2, double multiplier)
00439 {
00440 double window = 100.;
00441 if (which == 0) window = 2.*30.;
00442 else if (which == 1) window = 2.*30.;
00443 else if (which == 2) window = 2.*20.;
00444 else if (which == 3) window = 2.*50.;
00445
00446 TH1F *hist = dir->make<TH1F>(name.c_str(), "", 200, -window, window);
00447 for (std::vector<double*>::const_iterator r = residuals_begin(); r != residuals_end(); ++r)
00448 {
00449 double weight = 1./(*r)[whichredchi2];
00450 if (TMath::Prob(1./weight*12, 12) < 0.99) hist->Fill(multiplier * (*r)[which], weight);
00451 }
00452 }
00453
00454
00455 void MuonResidualsFitter::computeHistogramRangeAndBinning(int which, int &nbins, double &a, double &b)
00456 {
00457
00458 double *data = new double[numResiduals()];
00459 int n = 0;
00460 for (std::vector<double*>::const_iterator r = m_residuals.begin(); r != m_residuals.end(); r++)
00461 if (fabs((*r)[which])<50.)
00462 {
00463 data[n] = (*r)[which];
00464 n++;
00465 }
00466
00467
00468 const int n_quantiles = 7;
00469 double probabilities[n_quantiles] = {0.00135, 0.02275, 0.25, 0.5, 0.75, 0.97725, 0.99865};
00470
00471 double quantiles[n_quantiles];
00472 std::sort(data, data + n);
00473 TMath::Quantiles(n, n_quantiles, data, quantiles, probabilities, true, NULL, 7);
00474 delete [] data;
00475 double iqr = quantiles[4] - quantiles[2];
00476
00477
00478 double hbin = 2 * iqr / pow( n, 1./3);
00479
00480 a = quantiles[1];
00481 b = quantiles[5];
00482 nbins = (int) ( (b - a) / hbin + 3. );
00483
00484 std::cout<<" quantiles: "; for (int i=0;i<n_quantiles;i++) std::cout<<quantiles[i]<<" "; std::cout<<std::endl;
00485
00486
00487 std::cout<<" optimal h="<<hbin<<" nbins="<<nbins<<std::endl;
00488 }
00489
00490
00491 void MuonResidualsFitter::histogramChi2GaussianFit(int which, double &fit_mean, double &fit_sigma)
00492 {
00493 int nbins;
00494 double a, b;
00495 computeHistogramRangeAndBinning(which, nbins, a, b);
00496 if (a==b || a > b) { fit_mean = a; fit_sigma = 0; return; }
00497
00498 TH1D *hist = new TH1D("htmp", "", nbins, a, b);
00499 for (std::vector<double*>::const_iterator r = m_residuals.begin(); r != m_residuals.end(); ++r) hist->Fill( (*r)[which] );
00500
00501
00502 TF1 *f1= new TF1("f1","gaus", a, b);
00503 f1->SetParameter(0, hist->GetEntries());
00504 f1->SetParameter(1, 0);
00505 f1->SetParameter(2, hist->GetRMS());
00506 hist->Fit("f1","RQ");
00507
00508 fit_mean = f1->GetParameter(1);
00509 fit_sigma = f1->GetParameter(2);
00510 std::cout<<" h("<<nbins<<","<<a<<","<<b<<") mu="<<fit_mean<<" sig="<<fit_sigma<<std::endl;
00511
00512 delete f1;
00513 delete hist;
00514 }
00515
00516
00517
00518 void MuonResidualsFitter::selectPeakResiduals_simple(double nsigma, int nvar, int *vars)
00519 {
00520
00521 if (numResiduals()<25) return;
00522
00523 int nbefore = numResiduals();
00524
00525
00526 assert(nvar<=10);
00527
00528
00529 for (int v = 0; v<nvar; v++)
00530 {
00531 int which = vars[v];
00532 histogramChi2GaussianFit(which, m_center[which], m_radii[which]);
00533 m_radii[which] = nsigma * m_radii[which];
00534 }
00535
00536
00537 std::vector<double*>::iterator r = m_residuals.begin();
00538 while (r != m_residuals.end())
00539 {
00540 double ellipsoid_sum = 0;
00541 for (int v = 0; v<nvar; v++)
00542 {
00543 int which = vars[v];
00544 if (m_radii[which] == 0.) continue;
00545 ellipsoid_sum += pow( ( (*r)[which] - m_center[which]) / m_radii[which] , 2);
00546 }
00547 if (ellipsoid_sum <= 1.) ++r;
00548 else
00549 {
00550 delete [] (*r);
00551 r = m_residuals.erase(r);
00552 }
00553 }
00554 std::cout<<" N residuals "<<nbefore<<" -> "<<numResiduals()<<std::endl;
00555 }
00556
00557
00558
00559 void MuonResidualsFitter::selectPeakResiduals(double nsigma, int nvar, int *vars)
00560 {
00561
00562 for (int i=0; i<nvar; ++i) std::cout<<vars[i]<<" ";
00563 std::cout<<std::endl;
00564
00565
00566 if (numResiduals()<50) return;
00567
00568 size_t nbefore = numResiduals();
00569 std::cout<<" N residuals "<<nbefore<<" ~ "<<(size_t) std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true)<<std::endl;
00570
00571 assert(nvar<=10 && nvar>0);
00572
00573 std::vector<double*>::iterator r = m_residuals.begin();
00574
00575
00576 if (nvar==1)
00577 {
00578
00579 double *data = new double[nbefore];
00580 for (size_t i = 0; i < nbefore; i++) data[i] = m_residuals[i][ vars[0] ];
00581 double peak, sigma;
00582 TRobustEstimator re;
00583 re.EvaluateUni(nbefore, data, peak, sigma);
00584
00585
00586 while (r != m_residuals.end())
00587 {
00588 double distance = fabs( ((*r)[ vars[0] ] - peak)/sigma );
00589 if (distance <= nsigma) ++r;
00590 else
00591 {
00592 m_residuals_ok[r - m_residuals.begin()] = false;
00593 ++r;
00594
00595
00596 }
00597 }
00598 std::cout<<" N residuals "<<nbefore<<" -> "<<numResiduals()<<std::endl;
00599 return;
00600 }
00601
00602
00603 TRobustEstimator re(nbefore, nvar);
00604 r = m_residuals.begin();
00605 while (r != m_residuals.end())
00606 {
00607 double *row = new double[nvar];
00608 for (int v = 0; v<nvar; v++) row[v] = (*r)[ vars[v] ];
00609 re.AddRow(row);
00610 delete[] row;
00611 ++r;
00612 }
00613 re.Evaluate();
00614
00615
00616 TVectorD M(nvar);
00617 re.GetMean(M);
00618 TMatrixDSym Cov(nvar);
00619 re.GetCovariance(Cov);
00620 Cov.Invert();
00621
00622
00623 double conf_1d = TMath::Erf(nsigma/sqrt(2));
00624 double surf_radius = sqrt(TMath::ChisquareQuantile(conf_1d, nvar));
00625
00626
00627 r = m_residuals.begin();
00628 while (r != m_residuals.end())
00629 {
00630 TVectorD res(nvar);
00631 for (int v = 0; v<nvar; v++) res[v] = (*r)[ vars[v] ];
00632 double distance = sqrt( Cov.Similarity(res - M) );
00633 if (distance <= surf_radius) ++r;
00634 else
00635 {
00636 m_residuals_ok[r - m_residuals.begin()] = false;
00637 ++r;
00638
00639
00640 }
00641 }
00642 std::cout<<" N residuals "<<nbefore<<" -> "<<(size_t) std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true)<<std::endl;
00643 }
00644
00645
00646 void MuonResidualsFitter::correctBField(int idx_momentum, int idx_q)
00647 {
00648 const int Nbin = 17;
00649
00650 double min_pt = 9999.;
00651 for (std::vector<double*>::const_iterator r = residuals_begin(); r != residuals_end(); ++r)
00652 {
00653 double pt = fabs((*r)[idx_momentum]);
00654 if (pt < min_pt) min_pt = pt;
00655 }
00656 min_pt -= 0.01;
00657 const double bin_width = 1./min_pt/Nbin;
00658
00659
00660 std::vector<size_t> pos[Nbin], neg[Nbin], to_erase;
00661 for (size_t i = 0; i < m_residuals.size(); i++)
00662 {
00663 if (!m_residuals_ok[i]) continue;
00664 int bin = (int)floor(1./fabs(m_residuals[i][idx_momentum])/bin_width);
00665 if (m_residuals[i][idx_q] > 0) pos[bin].push_back(i);
00666 else neg[bin].push_back(i);
00667 }
00668
00669
00670 for (int j = 0; j < Nbin; j++)
00671 {
00672 size_t psize = pos[j].size();
00673 size_t nsize = neg[j].size();
00674 if (psize == nsize) continue;
00675
00676 std::set<int> idx_set;
00677 if (psize > nsize)
00678 {
00679 while (idx_set.size() < psize - nsize) idx_set.insert( gRandom->Integer(psize) );
00680 for (std::set<int>::iterator it = idx_set.begin() ; it != idx_set.end(); it++ ) to_erase.push_back(pos[j][*it]);
00681 }
00682 else
00683 {
00684 while (idx_set.size() < nsize - psize) idx_set.insert( gRandom->Integer(nsize) );
00685 for (std::set<int>::iterator it = idx_set.begin() ; it != idx_set.end(); it++ ) to_erase.push_back(neg[j][*it]);
00686 }
00687 }
00688
00689 std::sort(to_erase.begin(), to_erase.end(), std::greater<int>() );
00690 for (std::vector<size_t>::const_iterator e = to_erase.begin(); e != to_erase.end(); ++e)
00691 {
00692 m_residuals_ok[*e] = false;
00693
00694
00695 }
00696
00697 std::vector<size_t> apos[Nbin], aneg[Nbin];
00698 for (size_t i = 0; i < m_residuals.size(); i++)
00699 {
00700 if (!m_residuals_ok[i]) continue;
00701 int bin = (int)floor(1./fabs(m_residuals[i][idx_momentum])/bin_width);
00702 if (m_residuals[i][idx_q] > 0) apos[bin].push_back(i);
00703 else aneg[bin].push_back(i);
00704 }
00705 for (int j = 0; j < Nbin; j++) std::cout << "bin " << j << ": [pos,neg] sizes before & after: ["
00706 << pos[j].size() <<","<< neg[j].size()<<"] -> [" << apos[j].size() <<","<< aneg[j].size() << "]" << std::endl;
00707 std::cout<<" N residuals "<<m_residuals.size()<<" -> "<<(size_t) std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true)<<std::endl;
00708 }
00709
00710
00711 void MuonResidualsFitter::eraseNotSelectedResiduals()
00712 {
00713
00714 size_t n_ok = (size_t) std::count(m_residuals_ok.begin(), m_residuals_ok.end(), true);
00715 std::vector<double*> tmp(n_ok, 0);
00716 std::cout << "residuals sizes: all=" << m_residuals.size()<<" good="<<n_ok<<std::endl;
00717 int iok=0;
00718 for (size_t i = 0; i < m_residuals.size(); i++)
00719 {
00720 if (!m_residuals_ok[i]) continue;
00721 tmp[iok++] = m_residuals[i];
00722 }
00723 m_residuals.swap(tmp);
00724
00725 std::vector<bool> tmp_ok(n_ok, true);
00726 m_residuals_ok.swap(tmp_ok);
00727
00728 std::cout << "residuals size after eraseNotSelectedResiduals =" << m_residuals.size()<<" ok size="<<m_residuals_ok.size()<<std::endl;
00729 }