00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <TH1D.h>
00011 #include <TF1.h>
00012 #include <TGraph.h>
00013 #include <TString.h>
00014 #include <TPad.h>
00015 #include <TRandom3.h>
00016 #include <TLegend.h>
00017 #include <TMath.h>
00018 #include <Math/ProbFuncMathCore.h>
00019 #include <TString.h>
00020
00021 #include "../../BAT/BCLog.h"
00022 #include "../../BAT/BCDataSet.h"
00023 #include "../../BAT/BCDataPoint.h"
00024 #include "../../BAT/BCMath.h"
00025
00026 #include "BCHistogramFitter.h"
00027
00028
00029
00030 BCHistogramFitter::BCHistogramFitter()
00031 : BCModel()
00032 , fHistogram(0)
00033 , fFitFunction(0)
00034 , fHistogramExpected(0)
00035 {
00036
00037 MCMCSetNIterationsRun(2000);
00038 SetFillErrorBand();
00039 fFlagIntegration = true;
00040 flag_discrete = true;
00041 }
00042
00043
00044
00045 BCHistogramFitter::BCHistogramFitter(const char * name)
00046 : BCModel(name)
00047 , fHistogram(0)
00048 , fFitFunction(0)
00049 , fHistogramExpected(0)
00050 {
00051
00052 MCMCSetNIterationsRun(2000);
00053 SetFillErrorBand();
00054 fFlagIntegration = true;
00055 flag_discrete = true;
00056 }
00057
00058
00059
00060 BCHistogramFitter::BCHistogramFitter(TH1D * hist, TF1 * func)
00061 : BCModel()
00062 , fHistogram(0)
00063 , fFitFunction(0)
00064 , fHistogramExpected(0)
00065 {
00066 SetHistogram(hist);
00067 SetFitFunction(func);
00068
00069 MCMCSetNIterationsRun(2000);
00070 SetFillErrorBand();
00071 fFlagIntegration = true;
00072 flag_discrete = true;
00073 }
00074
00075
00076
00077 BCHistogramFitter::BCHistogramFitter(const char * name, TH1D * hist, TF1 * func)
00078 : BCModel(name)
00079 , fHistogram(0)
00080 , fFitFunction(0)
00081 , fHistogramExpected(0)
00082 {
00083 SetHistogram(hist);
00084 SetFitFunction(func);
00085
00086 MCMCSetNIterationsRun(2000);
00087 SetFillErrorBand();
00088 fFlagIntegration = true;
00089 flag_discrete = true;
00090 }
00091
00092
00093
00094 int BCHistogramFitter::SetHistogram(TH1D * hist)
00095 {
00096
00097 if (!hist) {
00098 BCLog::OutError("BCHistogramFitter::SetHistogram : TH1D not created.");
00099 return 0;
00100 }
00101
00102
00103 fHistogram = hist;
00104
00105
00106
00107
00108 BCDataSet * ds = new BCDataSet();
00109
00110
00111
00112
00113 int nbins = fHistogram->GetNbinsX();
00114 for (int i = 0; i < nbins; ++i) {
00115 BCDataPoint* dp = new BCDataPoint(2);
00116 dp ->SetValue(0, fHistogram->GetBinLowEdge(i + 1));
00117 dp ->SetValue(1, fHistogram->GetBinContent(i + 1));
00118 ds->AddDataPoint(dp);
00119 }
00120
00121
00122 SetDataSet(ds);
00123
00124
00125 double xmin = hist->GetBinLowEdge(1);
00126 double xmax = hist->GetBinLowEdge(nbins + 1);
00127
00128
00129 double histymin = hist->GetMinimum();
00130 double histymax = hist->GetMaximum();
00131
00132
00133
00134 double ymin = TMath::Max(0., histymin - 5. * sqrt(histymin));
00135 double ymax = histymax + 5. * sqrt(histymax);
00136
00137
00138 SetDataBoundaries(0, xmin, xmax);
00139 SetDataBoundaries(1, ymin, ymax);
00140
00141
00142 SetFitFunctionIndices(0, 1);
00143
00144
00145 return 1;
00146 }
00147
00148
00149
00150 int BCHistogramFitter::SetHistogramExpected(const std::vector<double> & parameters)
00151 {
00152
00153 if (fHistogramExpected) {
00154 delete fHistogramExpected;
00155 }
00156
00157 fHistogramExpected = new TH1D(*fHistogram);
00158
00159
00160 int nBins = fHistogramExpected->GetNbinsX();
00161
00162
00163 double dx = fHistogramExpected->GetBinWidth(1);
00164
00165
00166 fFitFunction->SetParameters(¶meters[0]);
00167
00168
00169 double fedgelow = fFitFunction->Eval(fHistogramExpected->GetBinLowEdge(1));
00170
00171
00172 for (int ibin = 1; ibin <= nBins; ++ibin) {
00173
00174 double xedgehi = fHistogramExpected->GetBinLowEdge(ibin) + dx;
00175
00176
00177 double yexp = 0.;
00178
00179
00180 if (fFlagIntegration)
00181 yexp = fFitFunction->Integral(xedgehi - dx, xedgehi);
00182
00183
00184 else {
00185
00186 double fedgehi = fFitFunction->Eval(xedgehi);
00187 yexp = fedgelow * dx + (fedgehi - fedgelow) * dx / 2.;
00188
00189
00190 fedgelow = fedgehi;
00191 }
00192
00193
00194 fHistogramExpected->SetBinContent(ibin, yexp);
00195
00196
00197 fHistogramExpected->SetBinError(ibin, 0.0);
00198
00199
00200 fHistogram->SetBinError(ibin, sqrt(yexp));
00201
00202 }
00203 return 1;
00204 }
00205
00206
00207
00208 int BCHistogramFitter::SetFitFunction(TF1 * func)
00209 {
00210
00211 if (!func) {
00212 BCLog::OutError("BCHistogramFitter::SetFitFunction : TF1 not created.");
00213 return 0;
00214 }
00215
00216
00217 fFitFunction = func;
00218
00219
00220 if(fName=="model")
00221 SetName(TString::Format("HistogramFitter with %s", fFitFunction->GetName()));
00222
00223
00224 fParameterSet->clear();
00225
00226
00227 int n = func->GetNpar();
00228
00229
00230 for (int i = 0; i < n; ++i) {
00231 double xmin;
00232 double xmax;
00233
00234 func->GetParLimits(i, xmin, xmax);
00235
00236 AddParameter(func->GetParName(i), xmin, xmax);
00237 }
00238
00239
00240 SetPriorConstantAll();
00241
00242 return GetNParameters();
00243 }
00244
00245
00246
00247 BCHistogramFitter::~BCHistogramFitter()
00248 {
00249 delete fHistogramExpected;
00250 }
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268 double BCHistogramFitter::LogLikelihood(const std::vector<double> & params)
00269 {
00270
00271 double loglikelihood = 0;
00272
00273
00274 fFitFunction->SetParameters(¶ms[0]);
00275
00276
00277 int nbins = fHistogram->GetNbinsX();
00278
00279
00280 double dx = fHistogram->GetBinWidth(1);
00281
00282
00283 double fedgelow = fFitFunction->Eval(fHistogram->GetBinLowEdge(1));
00284
00285
00286 for (int ibin = 1; ibin <= nbins; ++ibin) {
00287
00288 double xedgehi = fHistogram->GetBinLowEdge(ibin) + dx;
00289
00290
00291 double fedgehi = fFitFunction->Eval(xedgehi);
00292
00293
00294 double y = fHistogram->GetBinContent(ibin);
00295
00296 double yexp = 0.;
00297
00298
00299 if (fFlagIntegration)
00300 yexp = fFitFunction->Integral(xedgehi - dx, xedgehi);
00301
00302
00303 else {
00304 yexp = fedgelow * dx + (fedgehi - fedgelow) * dx / 2.;
00305
00306
00307 fedgelow = fedgehi;
00308 }
00309
00310
00311 loglikelihood += BCMath::LogPoisson(y, yexp);
00312 }
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327 return loglikelihood;
00328 }
00329
00330
00331
00332 double BCHistogramFitter::FitFunction(const std::vector<double> & x, const std::vector<double> & params)
00333 {
00334
00335 fFitFunction->SetParameters(¶ms[0]);
00336
00337 return fFitFunction->Eval(x[0]) * fHistogram->GetBinWidth(fHistogram->FindBin(x[0]));
00338 }
00339
00340
00341
00342 int BCHistogramFitter::Fit(TH1D * hist, TF1 * func)
00343 {
00344
00345 if (hist)
00346 SetHistogram(hist);
00347 else {
00348 BCLog::OutError("BCHistogramFitter::Fit : Histogram not defined.");
00349 return 0;
00350 }
00351
00352
00353 if (func)
00354 SetFitFunction(func);
00355 else {
00356 BCLog::OutError("BCHistogramFitter::Fit : Fit function not defined.");
00357 return 0;
00358 }
00359
00360 return Fit();
00361 }
00362
00363
00364
00365 int BCHistogramFitter::Fit()
00366 {
00367
00368 if (!fHistogram) {
00369 BCLog::OutError("BCHistogramFitter::Fit : Histogram not defined.");
00370 return 0;
00371 }
00372
00373
00374 if (!fFitFunction) {
00375 BCLog::OutError("BCHistogramFitter::Fit : Fit function not defined.");
00376 return 0;
00377 }
00378
00379
00380 MarginalizeAll();
00381
00382
00383
00384 FindModeMinuit(GetBestFitParameters(), -1);
00385
00386
00387 double pvalue;
00388 if (CalculatePValueFast(GetBestFitParameters(), pvalue))
00389 fPValue = pvalue;
00390 else
00391 BCLog::OutWarning(
00392 "BCHistogramFitter::Fit : Could not use the fast p-value evaluation.");
00393
00394
00395 PrintShortFitSummary();
00396
00397
00398 return 1;
00399 }
00400
00401
00402
00403 void BCHistogramFitter::DrawFit(const char * options, bool flaglegend)
00404 {
00405 if (!fHistogram) {
00406 BCLog::OutError("BCHistogramFitter::DrawFit : Histogram not defined.");
00407 return;
00408 }
00409
00410 if (!fFitFunction) {
00411 BCLog::OutError("BCHistogramFitter::DrawFit : Fit function not defined.");
00412 return;
00413 }
00414
00415 if (!fErrorBandXY || fBestFitParameters.empty()) {
00416 BCLog::OutError("BCHistogramFitter::DrawFit : Fit not performed yet.");
00417 return;
00418 }
00419
00420
00421 TString opt = options;
00422 opt.ToLower();
00423
00424
00425 if (!opt.Contains("same"))
00426 fHistogram->Draw(opt.Data());
00427
00428
00429 fErrorBand = GetErrorBandGraph(0.16, 0.84);
00430 fErrorBand->Draw("f same");
00431
00432
00433 fHistogram->Draw(TString::Format("%ssame", opt.Data()));
00434
00435
00436 fGraphFitFunction = GetFitFunctionGraph(GetBestFitParameters());
00437 fGraphFitFunction->SetLineColor(kRed);
00438 fGraphFitFunction->SetLineWidth(2);
00439 fGraphFitFunction->Draw("l same");
00440
00441
00442 if (flaglegend) {
00443 TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95);
00444 legend->SetBorderSize(0);
00445 legend->SetFillColor(kWhite);
00446 legend->AddEntry(fHistogram, "Data", "L");
00447 legend->AddEntry(fGraphFitFunction, "Best fit", "L");
00448 legend->AddEntry(fErrorBand, "Error band", "F");
00449 legend->Draw();
00450 }
00451
00452 gPad->RedrawAxis();
00453 }
00454
00455
00456 int BCHistogramFitter::CalculatePValueFast(const std::vector<double> & par, double &pvalue, int nIterations)
00457 {
00458
00459 return CalculatePValueFast(par, 0, pvalue, nIterations);
00460 }
00461
00462
00463 int BCHistogramFitter::CalculatePValueFast(const std::vector<double> & par,
00464 BCHistogramFitterToyDataInterface* callback, double &pvalue, int nIterations)
00465 {
00466
00467 if (par.size() != GetNParameters()) {
00468 BCLog::OutError(
00469 "BCHistogramFitter::CalculatePValueFast : Number of parameters is inconsistent.");
00470 return 0;
00471 }
00472
00473
00474 if (!fHistogram) {
00475 BCLog::OutError(
00476 "BCHistogramFitter::CalculatePValueFast : Histogram not defined.");
00477 return 0;
00478 }
00479
00480
00481 int nbins = fHistogram->GetNbinsX();
00482
00483 std::vector<int> histogram;
00484 std::vector<double> expectation;
00485 histogram.assign(nbins, 0);
00486 expectation.assign(nbins, 0);
00487
00488 double logp = 0;
00489 double logp_start = 0;
00490 int counter_pvalue = 0;
00491
00492
00493 SetHistogramExpected(par);
00494
00495
00496 for (int ibin = 0; ibin < nbins; ++ibin) {
00497
00498
00499 double yexp = fHistogramExpected->GetBinContent(ibin + 1);
00500
00501
00502 histogram[ibin] = int(yexp);
00503 expectation[ibin] = yexp;
00504
00505
00506 logp += BCMath::LogPoisson(double(int(yexp)), yexp);
00507
00508 logp_start += BCMath::LogPoisson(fHistogram->GetBinContent(ibin + 1),
00509 yexp);
00510 }
00511
00512
00513 for (int iiter = 0; iiter < nIterations; ++iiter) {
00514
00515 for (int ibin = 0; ibin < nbins; ++ibin) {
00516
00517 double ptest = fRandom->Rndm() - 0.5;
00518
00519
00520 if (ptest > 0) {
00521
00522 double r = expectation.at(ibin) / double(histogram.at(ibin) + 1);
00523
00524
00525 if (fRandom->Rndm() < r) {
00526 histogram[ibin] = histogram.at(ibin) + 1;
00527 logp += log(r);
00528 }
00529 }
00530
00531
00532 else if (ptest <= 0 && histogram[ibin] > 0) {
00533
00534 double r = double(histogram.at(ibin)) / expectation.at(ibin);
00535
00536
00537 if (fRandom->Rndm() < r) {
00538 histogram[ibin] = histogram.at(ibin) - 1;
00539 logp += log(r);
00540 }
00541 }
00542 }
00543
00544
00545
00546 if (callback)
00547 (*callback)(expectation, histogram);
00548
00549
00550 if (logp < logp_start)
00551 counter_pvalue++;
00552
00553 }
00554
00555
00556 pvalue = double(counter_pvalue) / double(nIterations);
00557
00558
00559 return 1;
00560 }
00561
00562
00563
00564 int BCHistogramFitter::CalculatePValueLikelihood(const std::vector<double> & par, double &pvalue)
00565 {
00566
00567 double logLambda = 0.0;
00568
00569
00570 SetHistogramExpected(par);
00571
00572 for (int ibin = 1; ibin <= fHistogram->GetNbinsX(); ++ibin) {
00573
00574
00575 double y = fHistogram->GetBinContent(ibin);
00576
00577
00578 double yexp = fHistogramExpected->GetBinContent(ibin);
00579
00580
00581 if (y == 0)
00582 logLambda += yexp;
00583 else
00584 logLambda += yexp - y + y * log(y / yexp);
00585 }
00586
00587
00588 logLambda *= 2.0;
00589
00590
00591 int nDoF = GetNDataPoints() - GetNParameters();
00592
00593
00594 fPValueNDoF = TMath::Prob(logLambda, nDoF);
00595 pvalue = fPValueNDoF;
00596
00597
00598 return 1;
00599 }
00600
00601
00602
00603 int BCHistogramFitter::CalculatePValueLeastSquares(const std::vector<double> & par, double &pvalue, bool weightExpect)
00604 {
00605
00606 double chi2 = 0.0;
00607
00608
00609 SetHistogramExpected(par);
00610
00611 for (int ibin = 1; ibin <= fHistogram->GetNbinsX(); ++ibin) {
00612
00613
00614 double y = fHistogram->GetBinContent(ibin);
00615
00616
00617 double yexp = fHistogramExpected->GetBinContent(ibin);
00618
00619
00620 double weight;
00621 if (weightExpect)
00622 weight = (yexp > 0) ? yexp : 1.0;
00623 else
00624 weight = (y > 0) ? y : 1.0;
00625
00626
00627 chi2 += (y - yexp) * (y - yexp) / weight;
00628 }
00629
00630
00631 int nDoF = GetNDataPoints() - GetNParameters();
00632
00633
00634 fPValueNDoF = TMath::Prob(chi2, nDoF);
00635 pvalue = fPValueNDoF;
00636
00637
00638 return 1;
00639 }
00640
00641
00642
00643 int BCHistogramFitter::CalculatePValueKolmogorov(const std::vector<double> & par, double& pvalue)
00644 {
00645 if (!fHistogramExpected || !fHistogram) {
00646 BCLog::OutError("BCHistogramFitter::CalculatePValueKolmogorov: "
00647 "Please define the reference distribution by calling \n"
00648 "BCHistogramFitter::SetHistogramExpected() first!");
00649 return 0;
00650 }
00651
00652
00653 SetHistogramExpected(par);
00654
00655
00656 pvalue = fHistogramExpected->KolmogorovTest(fHistogram);
00657
00658 fPValue = pvalue;
00659
00660
00661 return 1;
00662 }
00663
00664
00665
00666 double BCHistogramFitter::CDF(const std::vector<double>& parameters, int index, bool lower)
00667 {
00668
00669 if (parameters.empty())
00670 BCLog::OutWarning("BCHistogramFitter::CDF: parameter vector empty!");
00671
00672 index++;
00673
00674
00675 double yObs = fHistogram->GetBinContent(index);
00676
00677
00678
00679
00680
00681
00682 double edgeLow = fHistogram->GetBinLowEdge(index);
00683 double edgeHigh = edgeLow + fHistogram->GetBinWidth(index);
00684
00685
00686 double yExp = 0.0;
00687
00688
00689 if (fFlagIntegration)
00690 yExp = fFitFunction->Integral(edgeLow, edgeHigh, ¶meters[0]);
00691
00692
00693 else {
00694 double dx = fHistogram->GetBinWidth(index);
00695 double fEdgeLow = fFitFunction->Eval(edgeLow);
00696 double fEdgeHigh = fFitFunction->Eval(edgeHigh);
00697 yExp = fEdgeLow * dx + (fEdgeHigh - fEdgeLow) * dx / 2.;
00698 }
00699
00700
00701
00702
00703 if (lower) {
00704 if ((int) yObs >= 1)
00705 return ROOT::Math::poisson_cdf((unsigned int) (yObs - 1), yExp);
00706 else
00707 return 0.0;
00708 }
00709
00710 if ((double) (unsigned int) yObs != yObs) {
00711 BCLog::OutWarning(Form(
00712 "BCHistogramFitter::CDF: Histogram values should be integer!\n"
00713 " Bin %d = %f", index, yObs));
00714
00715
00716
00717
00718
00719 double U = fRandom->Rndm();
00720 double yObsLower = (unsigned int) yObs;
00721 if (U > (yObs - yObsLower))
00722 yObs = yObsLower;
00723 else
00724 yObs = yObsLower + 1;
00725 }
00726
00727
00728
00729
00730 return ROOT::Math::poisson_cdf((unsigned int) yObs, yExp);
00731 }
00732
00733