ROOT
6.07/01
Reference Guide
|
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 polynomails 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 emperical 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 24 of file BernsteinCorrection.h.
Public Member Functions | |
BernsteinCorrection (double tolerance=0.05) | |
virtual | ~BernsteinCorrection () |
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) |
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... | |
Private Attributes | |
Int_t | fMaxDegree |
Double_t | fMaxCorrection |
Double_t | fTolerance |
#include <RooStats/BernsteinCorrection.h>
BernsteinCorrection::BernsteinCorrection | ( | double | tolerance = 0.05 | ) |
Definition at line 87 of file BernsteinCorrection.cxx.
|
inlinevirtual |
Definition at line 28 of file BernsteinCorrection.h.
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 212 of file BernsteinCorrection.cxx.
Referenced by rs_bernsteinCorrection().
Int_t BernsteinCorrection::ImportCorrectedPdf | ( | RooWorkspace * | wks, |
const char * | nominalName, | ||
const char * | varName, | ||
const char * | dataName | ||
) |
Main method for Bernstein correction.
Definition at line 95 of file BernsteinCorrection.cxx.
Referenced by rs_bernsteinCorrection().
Definition at line 31 of file BernsteinCorrection.h.
Definition at line 32 of file BernsteinCorrection.h.
|
private |
Definition at line 44 of file BernsteinCorrection.h.
Referenced by CreateQSamplingDist(), ImportCorrectedPdf(), and SetMaxCorrection().
|
private |
Definition at line 43 of file BernsteinCorrection.h.
Referenced by ImportCorrectedPdf(), and SetMaxDegree().
|
private |
Definition at line 45 of file BernsteinCorrection.h.
Referenced by ImportCorrectedPdf().