The idea is that one has a distribution coming either from data or Monte Carlo (called "reality" in the macro) and a nominal model that is not sufficiently flexible to take into account the real distribution. One wants to take into account the systematic associated with this imperfect modeling by augmenting the nominal model with some correction term (in this case a polynomial). The BernsteinCorrection utility will import into your workspace a corrected model given by nominal(x) * poly_N(x), where poly_N is an n-th order polynomial in the Bernstein basis. The degree N of the polynomial is chosen by specifying the tolerance one has in adding an extra term to the polynomial. The Bernstein basis is nice because it only has positive-definite terms and works well with PDFs. Finally, the macro makes a plot of:
[#0] WARNING:InputArguments -- The parameter 'sigma' with range [0, 10] of the RooGaussian 'nominal' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing dataset realityData
[#1] INFO:ObjectHandling -- RooWorkSpace::import(myWorksspace) changing name of dataset from realityData to data
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::x
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooGaussian::nominal
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooConstVar::0
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::sigma
BernsteinCorrection::ImportCorrectedPdf - Doing initial Fit with nominal model
[#1] INFO:Fitting -- RooAbsPdf::fitTo(nominal_over_nominal_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_nominal_over_nominal_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_data) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooEffProd::corrected
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooBernstein::poly
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_0
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_1
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_2
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_3
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_4
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_5
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWorksspace) importing RooRealVar::c_6
------ Begin Bernstein Correction Log --------
degree = 1 -log L(0) = 1216.78 -log L(1) = 1208.89 q = 15.7692 P(chi^2_1 > q) = 7.1557e-05
degree = 2 -log L(1) = 1208.89 -log L(2) = 1203.21 q = 11.3692 P(chi^2_1 > q) = 0.000746732
degree = 3 -log L(2) = 1203.21 -log L(3) = 1198.85 q = 8.72171 P(chi^2_1 > q) = 0.00314444
degree = 4 -log L(3) = 1198.85 -log L(4) = 1190.16 q = 17.3777 P(chi^2_1 > q) = 3.06393e-05
degree = 5 -log L(4) = 1190.16 -log L(5) = 1183.56 q = 13.1965 P(chi^2_1 > q) = 0.00028048
degree = 6 -log L(5) = 1183.56 -log L(6) = 1182.57 q = 1.98352 P(chi^2_1 > q) = 0.15902
------ End Bernstein Correction Log --------
Correction based on Bernstein Poly of degree 6
[#1] INFO:Fitting -- RooAbsPdf::fitTo(nominal_over_nominal_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_nominal_over_nominal_Int[x]_realityData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
Minuit2Minimizer: Minimize with max-calls 500 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = 1216.77793416266286
Edm = 4.16187344562173991e-07
Nfcn = 19
sigma = 1.18138 +/- 0.0315451 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(corrected_over_corrected_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_corrected_over_corrected_Int[x]_realityData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
Minuit2Minimizer: Minimize with max-calls 3500 convergence for edm < 1 strategy 1
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
Minuit2Minimizer : Valid minimum - status = 0
FVAL = 1182.56856162134386
Edm = 0.000128904133136698251
Nfcn = 185
c_1 = 3.18336 +/- 0.783711 (limited)
c_2 = 4.80577e-06 +/- 3.11001 (limited)
c_3 = 1.19142e-06 +/- 1.51771 (limited)
c_4 = 0.966777 +/- 2.14189 (limited)
c_5 = 0.21563 +/- 67.5033 (limited)
c_6 = 10.4446 +/- 20.5651 (limited)
sigma = 1.26666 +/- 0.209609 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(corrected_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
made pdfs, make toy generator
on toy 0
on toy 1
on toy 2
on toy 3
on toy 4
on toy 5
on toy 6
on toy 7
on toy 8
on toy 9
on toy 10
on toy 11
on toy 12
on toy 13
on toy 14
on toy 15
on toy 16
on toy 17
on toy 18
on toy 19
#include <string>
#include <vector>
#include <stdio.h>
#include <sstream>
#include <iostream>
{
RooGaussian narrow(
"narrow",
"",
x, RooConst(0.), RooConst(.8));
RooGaussian wide(
"wide",
"",
x, RooConst(0.), RooConst(2.));
std::unique_ptr<RooDataSet>
data{reality.generate(
x, 1000)};
}
Int_t degree = bernsteinCorrection.ImportCorrectedPdf(wks,
"nominal",
"x",
"data");
if (degree < 0) {
Error(
"rs_bernsteinCorrection",
"Bernstein correction failed ! ");
return;
}
cout << " Correction based on Bernstein Poly of degree " << degree << endl;
nominal.plotOn(frame);
if (!corrected)
return;
if (poly)
bool checkSamplingDist = true;
int numToyMC = 20;
if (checkSamplingDist) {
}
if (checkSamplingDist) {
TH1F *samplingDist =
new TH1F(
"samplingDist",
"", 20, 0, 10);
TH1F *samplingDistExtra =
new TH1F(
"samplingDistExtra",
"", 20, 0, 10);
bernsteinCorrection.CreateQSamplingDist(wks, "nominal", "x", "data", samplingDist, samplingDistExtra, degree,
numToyMC);
samplingDistExtra->
Draw();
samplingDist->
Draw(
"same");
}
}
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
static void SetDefaultPrintLevel(int level)
Set the default Print Level.
Abstract interface for all probability density functions.
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Efficient implementation of a sum of PDFs of the form.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Plot frame and a container for graphics objects within that frame.
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Variable that can be changed from the outside.
BernsteinCorrection is a utility in RooStats to augment a nominal PDF with a polynomial correction te...
Persistable container for RooFit projects.
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
bool import(const RooAbsArg &arg, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={})
Import a RooAbsArg object, e.g.
virtual void SetLineColor(Color_t lcolor)
Set the line color.
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
1-D histogram with a float per channel (see TH1 documentation)
void Draw(Option_t *option="") override
Draw this histogram with options.
RooCmdArg PrintLevel(Int_t code)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.