00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <TH1D.h>
00011 #include <TH2D.h>
00012 #include <TF1.h>
00013 #include <TGraph.h>
00014 #include <TGraphAsymmErrors.h>
00015 #include <TString.h>
00016 #include <TPad.h>
00017 #include <TRandom3.h>
00018 #include <TLegend.h>
00019 #include <TMath.h>
00020
00021 #include "../../BAT/BCLog.h"
00022 #include "../../BAT/BCDataSet.h"
00023 #include "../../BAT/BCDataPoint.h"
00024 #include "../../BAT/BCMath.h"
00025 #include "../../BAT/BCH1D.h"
00026
00027 #include "BCEfficiencyFitter.h"
00028
00029
00030
00031 BCEfficiencyFitter::BCEfficiencyFitter()
00032 : BCModel()
00033 , fHistogram1(0)
00034 , fHistogram2(0)
00035 , fFitFunction(0)
00036 , fHistogramBinomial(0)
00037 , fDataPointType(1)
00038 {
00039
00040 MCMCSetNIterationsRun(2000);
00041 MCMCSetRValueCriterion(0.01);
00042 SetFillErrorBand();
00043 fFlagIntegration = false;
00044 }
00045
00046
00047
00048 BCEfficiencyFitter::BCEfficiencyFitter(const char * name)
00049 : BCModel(name)
00050 , fHistogram1(0)
00051 , fHistogram2(0)
00052 , fFitFunction(0)
00053 , fHistogramBinomial(0)
00054 , fDataPointType(1)
00055 {
00056
00057 MCMCSetNIterationsRun(2000);
00058 MCMCSetRValueCriterion(0.01);
00059 SetFillErrorBand();
00060 fFlagIntegration = false;
00061 }
00062
00063
00064
00065 BCEfficiencyFitter::BCEfficiencyFitter(TH1D * hist1, TH1D * hist2, TF1 * func)
00066 : BCModel()
00067 , fHistogram1(0)
00068 , fHistogram2(0)
00069 , fFitFunction(0)
00070 , fHistogramBinomial(0)
00071 , fDataPointType(1)
00072 {
00073 SetHistograms(hist1, hist2);
00074 SetFitFunction(func);
00075
00076 MCMCSetNIterationsRun(2000);
00077 MCMCSetRValueCriterion(0.01);
00078 SetFillErrorBand();
00079 fFlagIntegration = false;
00080 }
00081
00082
00083
00084 BCEfficiencyFitter::BCEfficiencyFitter(const char * name, TH1D * hist1, TH1D * hist2, TF1 * func)
00085 : BCModel(name)
00086 , fHistogram1(0)
00087 , fHistogram2(0)
00088 , fFitFunction(0)
00089 , fHistogramBinomial(0)
00090 , fDataPointType(1)
00091 {
00092 SetHistograms(hist1, hist2);
00093 SetFitFunction(func);
00094
00095 MCMCSetNIterationsRun(2000);
00096 MCMCSetRValueCriterion(0.01);
00097 SetFillErrorBand();
00098 fFlagIntegration = false;
00099 }
00100
00101
00102
00103 int BCEfficiencyFitter::SetHistograms(TH1D * hist1, TH1D * hist2)
00104 {
00105
00106 if (!hist1 || !hist2) {
00107 BCLog::OutError("BCEfficiencyFitter::SetHistograms : TH1D not created.");
00108 return 0;
00109 }
00110
00111
00112 if (hist1->GetNbinsX() != hist2->GetNbinsX()) {
00113 BCLog::OutError("BCEfficiencyFitter::SetHistograms : Histograms do not have the same number of bins.");
00114 return 0;
00115 }
00116
00117
00118 for (int i = 1; i <= hist1->GetNbinsX(); ++i)
00119 if (hist1->GetBinContent(i) < hist2->GetBinContent(i)) {
00120 BCLog::OutError("BCEfficiencyFitter::SetHistograms : Histogram 1 has fewer entries than histogram 2.");
00121 return 0;
00122 }
00123
00124
00125 fHistogram1 = hist1;
00126 fHistogram2 = hist2;
00127
00128
00129
00130
00131 BCDataSet * ds = new BCDataSet();
00132
00133
00134 int nbins = fHistogram1->GetNbinsX();
00135 for (int i = 0; i < nbins; ++i) {
00136 BCDataPoint* dp = new BCDataPoint(2);
00137 ds->AddDataPoint(dp);
00138 }
00139
00140
00141 SetDataSet(ds);
00142
00143
00144 double xmin = hist1->GetBinLowEdge(1);
00145 double xmax = hist1->GetBinLowEdge(nbins+1);
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157 SetDataBoundaries(0, xmin, xmax);
00158 SetDataBoundaries(1, 0.0, 1.0);
00159
00160
00161 SetFitFunctionIndices(0, 1);
00162
00163
00164 return 1;
00165 }
00166
00167
00168
00169 int BCEfficiencyFitter::SetFitFunction(TF1 * func)
00170 {
00171
00172 if(!func) {
00173 BCLog::OutError("BCEfficiencyFitter::SetFitFunction : TF1 not created.");
00174 return 0;
00175 }
00176
00177
00178 fFitFunction = func;
00179
00180
00181 if(fName=="model")
00182 SetName(TString::Format("BCEfficiencyFitter with %s",fFitFunction->GetName()));
00183
00184
00185 fParameterSet->clear();
00186
00187
00188 int n = func->GetNpar();
00189
00190
00191 for (int i = 0; i < n; ++i)
00192 {
00193 double xmin;
00194 double xmax;
00195
00196 func->GetParLimits(i, xmin, xmax);
00197
00198 AddParameter(func->GetParName(i), xmin, xmax);
00199 }
00200
00201
00202 SetPriorConstantAll();
00203
00204 return GetNParameters();;
00205 }
00206
00207
00208
00209 BCEfficiencyFitter::~BCEfficiencyFitter()
00210 {
00211 if (fHistogram1)
00212 delete fHistogram1;
00213
00214 if (fHistogram2)
00215 delete fHistogram2;
00216
00217 if (fHistogramBinomial)
00218 delete fHistogramBinomial;
00219 }
00220
00221
00222
00223 void BCEfficiencyFitter::SetDataPointType(int type)
00224 {
00225 if(type < 0 || type > 2)
00226 BCLog::OutError(Form("BCEfficiencyFitter::SetDataPointType : Unknown data point type %d (should be between 0 and 2).",type));
00227 else
00228 fDataPointType = type;
00229 }
00230
00231
00232
00233 double BCEfficiencyFitter::LogLikelihood(const std::vector<double> & params)
00234 {
00235
00236
00237 double loglikelihood = 0;
00238
00239
00240 fFitFunction->SetParameters(¶ms[0]);
00241
00242
00243 int nbins = fHistogram1->GetNbinsX();
00244
00245
00246 double dx = fHistogram1->GetXaxis()->GetBinWidth(0);
00247
00248
00249 for (int ibin = 1; ibin <= nbins; ++ibin)
00250 {
00251
00252 int n = int(fHistogram1->GetBinContent(ibin));
00253
00254
00255 int k = int(fHistogram2->GetBinContent(ibin));
00256
00257
00258 double x = fHistogram1->GetBinCenter(ibin);
00259
00260 double eff = 0;
00261
00262
00263 if (fFlagIntegration)
00264 eff = fFitFunction->Integral(x - dx/2., x + dx/2.) / dx;
00265
00266
00267 else
00268 eff = (fFitFunction->Eval(x + dx/2.) + fFitFunction->Eval(x - dx/2.)) / 2.;
00269
00270
00271 loglikelihood += BCMath::LogApproxBinomial(n, k, eff);
00272 }
00273
00274 return loglikelihood;
00275 }
00276
00277
00278
00279 double BCEfficiencyFitter::FitFunction(const std::vector<double> & x, const std::vector<double> & params)
00280 {
00281
00282 fFitFunction->SetParameters(¶ms[0]);
00283
00284 return fFitFunction->Eval(x[0]);
00285 }
00286
00287
00288
00289 int BCEfficiencyFitter::Fit(TH1D * hist1, TH1D * hist2, TF1 * func)
00290 {
00291
00292 if (hist1 && hist2)
00293 SetHistograms(hist1, hist2);
00294 else {
00295 BCLog::OutError("BCEfficiencyFitter::Fit : Histogram(s) not defined.");
00296 return 0;
00297 }
00298
00299
00300 if (func)
00301 SetFitFunction(func);
00302 else {
00303 BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
00304 return 0;
00305 }
00306
00307 return Fit();
00308 }
00309
00310
00311
00312 int BCEfficiencyFitter::Fit()
00313 {
00314
00315 if (!fHistogram1 || !fHistogram2) {
00316 BCLog::OutError("BCEfficiencyFitter::Fit : Histogram(s) not defined.");
00317 return 0;
00318 }
00319
00320
00321 if (!fFitFunction) {
00322 BCLog::OutError("BCEfficiencyFitter::Fit : Fit function not defined.");
00323 return 0;
00324 }
00325
00326
00327 MarginalizeAll();
00328
00329
00330
00331 FindModeMinuit( GetBestFitParameters() , -1);
00332
00333
00334 double pvalue;
00335 if ( CalculatePValueFast(GetBestFitParameters(), pvalue) )
00336 fPValue = pvalue;
00337 else
00338 BCLog::OutError("BCEfficiencyFitter::Fit : Could not use the fast p-value evaluation.");
00339
00340
00341 PrintShortFitSummary();
00342
00343
00344 return 1;
00345 }
00346
00347
00348
00349 void BCEfficiencyFitter::DrawFit(const char * options, bool flaglegend)
00350 {
00351 if (!fHistogram1 || !fHistogram2) {
00352 BCLog::OutError("BCEfficiencyFitter::DrawFit : Histogram(s) not defined.");
00353 return;
00354 }
00355
00356 if (!fFitFunction) {
00357 BCLog::OutError("BCEfficiencyFitter::DrawFit : Fit function not defined.");
00358 return;
00359 }
00360
00361
00362 TGraphAsymmErrors * histRatio = new TGraphAsymmErrors();
00363 histRatio->SetMarkerStyle(20);
00364 histRatio->SetMarkerSize(1.5);
00365
00366 int npoints = 0;
00367
00368
00369 for (int i = 1; i <= fHistogram1->GetNbinsX(); ++i)
00370 {
00371 if(int(fHistogram1->GetBinContent(i)) == 0) {
00372 ++npoints;
00373 continue;
00374 }
00375
00376
00377 double xexp, xmin, xmax;
00378 GetUncertainties( int(fHistogram1->GetBinContent(i)), int(fHistogram2->GetBinContent(i)), 0.68, xexp, xmin, xmax);
00379
00380 histRatio->SetPoint( npoints, fHistogram1->GetBinCenter(i), xexp);
00381
00382
00383 histRatio->SetPointEXhigh(npoints, 0.);
00384 histRatio->SetPointEXlow(npoints, 0.);
00385
00386
00387 histRatio->SetPointEYhigh(npoints, xmax - xexp);
00388 histRatio->SetPointEYlow(npoints, xexp - xmin);
00389
00390 ++npoints;
00391 }
00392
00393
00394
00395 TString opt = options;
00396 opt.ToLower();
00397
00398
00399 if(!opt.Contains("same"))
00400 {
00401
00402 TH2D * hist_axes = new TH2D("hist_axes",
00403 Form(";%s;ratio", fHistogram1->GetXaxis()->GetTitle()),
00404 fHistogram1->GetNbinsX(),
00405 fHistogram1->GetXaxis()->GetBinLowEdge(1),
00406 fHistogram1->GetXaxis()->GetBinLowEdge(fHistogram1->GetNbinsX()+1),
00407 1, 0., 1.);
00408 hist_axes->SetStats(kFALSE);
00409 hist_axes->Draw();
00410
00411 histRatio->Draw(TString::Format("%sp",opt.Data()));
00412 }
00413
00414
00415 fErrorBand = GetErrorBandGraph(0.16, 0.84);
00416 fErrorBand->Draw("f same");
00417
00418
00419 histRatio->SetMarkerSize(.7);
00420 histRatio->Draw(TString::Format("%ssamep",opt.Data()));
00421
00422
00423 fGraphFitFunction = GetFitFunctionGraph( GetBestFitParameters() );
00424 fGraphFitFunction->SetLineColor(kRed);
00425 fGraphFitFunction->SetLineWidth(2);
00426 fGraphFitFunction->Draw("l same");
00427
00428
00429 if (flaglegend)
00430 {
00431 TLegend * legend = new TLegend(0.25, 0.75, 0.55, 0.95);
00432 legend->SetBorderSize(0);
00433 legend->SetFillColor(kWhite);
00434 legend->AddEntry(histRatio, "Data", "L");
00435 legend->AddEntry(fGraphFitFunction, "Best fit", "L");
00436 legend->AddEntry(fErrorBand, "Error band", "F");
00437 legend->Draw();
00438 }
00439 gPad->RedrawAxis();
00440 }
00441
00442
00443 int BCEfficiencyFitter::CalculatePValueFast(const std::vector<double> & par, double &pvalue)
00444 {
00445
00446 if (par.size() != GetNParameters()) {
00447 BCLog::OutError("BCEfficiencyFitter::CalculatePValueFast : Number of parameters is inconsistent.");
00448 return 0;
00449 }
00450
00451
00452 if (!fHistogram1 || !fHistogram2) {
00453 BCLog::OutError("BCEfficiencyFitter::CalculatePValueFast : Histogram not defined.");
00454 return 0;
00455 }
00456
00457
00458 int nbins = fHistogram1->GetNbinsX();
00459
00460 std::vector<int> histogram;
00461 std::vector<double> expectation;
00462 histogram.assign(nbins, 0);
00463 expectation.assign(nbins, 0);
00464
00465 double logp = 0;
00466 double logp_start = 0;
00467 int counter_pvalue = 0;
00468
00469
00470 for (int ibin = 0; ibin < nbins; ++ibin) {
00471
00472 double xmin = fHistogram1->GetBinLowEdge(ibin+1);
00473 double xmax = fHistogram1->GetBinLowEdge(ibin+2);
00474
00475
00476 double yexp = fFitFunction->Integral(xmin, xmax);
00477
00478
00479 int n = int(fHistogram1->GetBinContent(ibin));
00480
00481
00482 int k = int(fHistogram2->GetBinContent(ibin));
00483
00484 histogram[ibin] = k;
00485 expectation[ibin] = n * yexp;
00486
00487
00488 logp += BCMath::LogApproxBinomial(n, k, yexp);
00489 }
00490 logp_start = logp;
00491
00492 int niter = 100000;
00493
00494
00495 for (int iiter = 0; iiter < niter; ++iiter)
00496 {
00497
00498 for (int ibin = 0; ibin < nbins; ++ibin)
00499 {
00500
00501 int n = int(fHistogram1->GetBinContent(ibin));
00502
00503
00504 int k = histogram.at(ibin);
00505
00506
00507 double yexp = 0;
00508 if (n > 0)
00509 yexp = expectation.at(ibin)/n;
00510
00511
00512 double ptest = fRandom->Rndm() - 0.5;
00513
00514
00515 if (!(yexp > 0. || yexp < 1.0))
00516 continue;
00517
00518
00519 if (ptest > 0 && (k < n))
00520 {
00521
00522 double r = (double(n)-double(k))/(double(k)+1.) * yexp / (1. - yexp);
00523
00524
00525 if (fRandom->Rndm() < r)
00526 {
00527 histogram[ibin] = k + 1;
00528 logp += log(r);
00529 }
00530 }
00531
00532
00533 else if (ptest <= 0 && (k > 0))
00534 {
00535
00536 double r = double(k) / (double(n)-(double(k)-1)) * (1-yexp)/yexp;
00537
00538
00539 if (fRandom->Rndm() < r)
00540 {
00541 histogram[ibin] = k - 1;
00542 logp += log(r);
00543 }
00544 }
00545
00546 }
00547
00548
00549 if (logp < logp_start)
00550 counter_pvalue++;
00551
00552 }
00553
00554
00555 pvalue = double(counter_pvalue) / double(niter);
00556
00557
00558 return 1;
00559 }
00560
00561
00562 int BCEfficiencyFitter::GetUncertainties(int n, int k, double p, double &xexp, double &xmin, double &xmax)
00563 {
00564 if (n == 0)
00565 {
00566 xexp = 0.;
00567 xmin = 0.;
00568 xmax = 0.;
00569 return 0;
00570 }
00571
00572 BCLog::OutDebug(Form("Calculating efficiency data-point of type %d for (n,k) = (%d,%d)",fDataPointType,n,k));
00573
00574
00575 if (fHistogramBinomial)
00576 fHistogramBinomial->Reset();
00577 else
00578 fHistogramBinomial = new TH1D("hist_binomial", "", 1000, 0., 1.);
00579
00580
00581 for (int i = 1; i <= 1000; ++i) {
00582 double x = fHistogramBinomial->GetBinCenter(i);
00583 double val = BCMath::ApproxBinomial(n, k, x);
00584 fHistogramBinomial->SetBinContent(i, val);
00585 }
00586
00587
00588 fHistogramBinomial->Scale(1. / fHistogramBinomial->Integral());
00589
00590
00591 if (fDataPointType == 0) {
00592 xexp = fHistogramBinomial->GetMean();
00593 double rms = fHistogramBinomial->GetRMS();
00594 xmin = xexp-rms;
00595 xmax = xexp+rms;
00596 BCLog::OutDebug(Form(" - mean = %f , rms = %f)",xexp,rms));
00597 }
00598 else if (fDataPointType == 1) {
00599 xexp = (double)k/(double)n;
00600 BCH1D * fbh = new BCH1D((TH1D*)fHistogramBinomial->Clone("hcp"));
00601 std::vector<double> intervals = fbh->GetSmallestIntervals(p);
00602 int ninter = intervals.size();
00603 if ( ninter<2 ) {
00604 xmin = xmax = xexp = 0.;
00605 }
00606 else {
00607 xmin = intervals[0];
00608 xmax = intervals[1];
00609 }
00610 delete fbh;
00611 }
00612 else {
00613
00614 int nprobSum = 3;
00615 double q[3];
00616 double probSum[3];
00617 probSum[0] = (1.-p)/2.;
00618 probSum[1] = 1.-(1.-p)/2.;
00619 probSum[2] = .5;
00620
00621 fHistogramBinomial->GetQuantiles(nprobSum, q, probSum);
00622
00623 xexp = q[2];
00624 xmin = q[0];
00625 xmax = q[1];
00626
00627 }
00628
00629 BCLog::OutDebug(Form(" - efficiency = %f , range (%f - %f)",xexp,xmin,xmax));
00630
00631 return 1;
00632 }
00633
00634