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