Logo ROOT  
Reference Guide
RooStats::BernsteinCorrection Class Reference

BernsteinCorrection is a utility in RooStats to augment a nominal PDF with a polynomial correction term.

This is useful for incorporating systematic variations to the nominal PDF. The Bernstein basis polynomials are particularly appropriate because they are positive definite.

This tool was inspired by the work of Glen Cowan together with Stephan Horner, Sascha Caron, Eilam Gross, and others. The initial implementation is independent work. The major step forward in the approach was to provide a well defined algorithm that specifies the order of polynomial to be included in the correction. This is an empirical algorithm, so in addition to the nominal model it needs either a real data set or a simulated one. In the early work, the nominal model was taken to be a histogram from Monte Carlo simulations, but in this implementation it is generalized to an arbitrary PDF (which includes a RooHistPdf). The algorithm basically consists of a hypothesis test of an nth-order correction (null) against a n+1-th order correction (alternate). The quantity q = -2 log LR is used to determine whether the n+1-th order correction is a major improvement to the n-th order correction. The distribution of q is expected to be roughly \(\chi^2\) with one degree of freedom if the n-th order correction is a good model for the data. Thus, one only moves to the n+1-th order correction of q is relatively large. The chance that one moves from the n-th to the n+1-th order correction when the n-th order correction (eg. a type 1 error) is sufficient is given by the Prob( \(\chi^2_1\) > threshold). The constructor of this class allows you to directly set this tolerance (in terms of probability that the n+1-th term is added unnecessarily).

To do: Add another method to the utility that will make the sampling distribution for -2 log lambda for various m vs. m+1 order corrections using a nominal model and perhaps having two ways of generating the toys (either via a histogram or via an independent model that is supposed to reflect reality). That will allow one to make plots like Glen has at the end of his DRAFT very easily.

Definition at line 22 of file BernsteinCorrection.h.

Public Member Functions

 BernsteinCorrection (double tolerance=0.05)
 
virtual ~BernsteinCorrection ()
 
void CreateQSamplingDist (RooWorkspace *wks, const char *nominalName, const char *varName, const char *dataName, TH1F *, TH1F *, Int_t degree, Int_t nToys=500)
 Create sampling distribution for q given degree-1 vs. degree corrections. More...
 
Int_t ImportCorrectedPdf (RooWorkspace *, const char *, const char *, const char *)
 Main method for Bernstein correction. More...
 
void SetMaxCorrection (Double_t maxCorr)
 
void SetMaxDegree (Int_t maxDegree)
 

Private Attributes

Double_t fMaxCorrection
 
Int_t fMaxDegree
 
Double_t fTolerance
 

#include <RooStats/BernsteinCorrection.h>

Constructor & Destructor Documentation

◆ BernsteinCorrection()

BernsteinCorrection::BernsteinCorrection ( double  tolerance = 0.05)

Definition at line 81 of file BernsteinCorrection.cxx.

◆ ~BernsteinCorrection()

virtual RooStats::BernsteinCorrection::~BernsteinCorrection ( )
inlinevirtual

Definition at line 26 of file BernsteinCorrection.h.

Member Function Documentation

◆ CreateQSamplingDist()

void BernsteinCorrection::CreateQSamplingDist ( RooWorkspace wks,
const char *  nominalName,
const char *  varName,
const char *  dataName,
TH1F samplingDist,
TH1F samplingDistExtra,
Int_t  degree,
Int_t  nToys = 500 
)

Create sampling distribution for q given degree-1 vs. degree corrections.

Definition at line 206 of file BernsteinCorrection.cxx.

◆ ImportCorrectedPdf()

Int_t BernsteinCorrection::ImportCorrectedPdf ( RooWorkspace wks,
const char *  nominalName,
const char *  varName,
const char *  dataName 
)

Main method for Bernstein correction.

Definition at line 89 of file BernsteinCorrection.cxx.

◆ SetMaxCorrection()

void RooStats::BernsteinCorrection::SetMaxCorrection ( Double_t  maxCorr)
inline

Definition at line 29 of file BernsteinCorrection.h.

◆ SetMaxDegree()

void RooStats::BernsteinCorrection::SetMaxDegree ( Int_t  maxDegree)
inline

Definition at line 30 of file BernsteinCorrection.h.

Member Data Documentation

◆ fMaxCorrection

Double_t RooStats::BernsteinCorrection::fMaxCorrection
private

Definition at line 42 of file BernsteinCorrection.h.

◆ fMaxDegree

Int_t RooStats::BernsteinCorrection::fMaxDegree
private

Definition at line 41 of file BernsteinCorrection.h.

◆ fTolerance

Double_t RooStats::BernsteinCorrection::fTolerance
private

Definition at line 43 of file BernsteinCorrection.h.


The documentation for this class was generated from the following files: