51 #ifndef RooStats_BernsteinCorrection 87 BernsteinCorrection::BernsteinCorrection(
Double_t tolerance):
88 fMaxDegree(10), fMaxCorrection(100), fTolerance(tolerance){
96 const char* nominalName,
98 const char* dataName){
104 if (!x || !nominal || !data) {
105 cout <<
"Error: wrong name for pdf or variable or dataset - return -1 " << std::endl;
109 std::cout <<
"BernsteinCorrection::ImportCorrectedPdf - Doing initial Fit with nominam model " << std::endl;
118 if (nominalResult->
status() != 0 ) {
119 std::cout <<
"BernsteinCorrection::ImportCorrectedPdf - Error fit with nominal model failed - exit" << std::endl;
124 std::stringstream
log;
125 log <<
"------ Begin Bernstein Correction Log --------" << endl;
129 vector<RooRealVar*> coefficients;
134 bool keepGoing =
true;
139 std::stringstream str;
143 "Bernstein basis poly coefficient",
146 coefficients.push_back(newCoef);
163 if (result->
status() != 0) {
164 std::cout <<
"BernsteinCorrection::ImportCorrectedPdf - Error fit with corrected model failed" << std::endl;
171 q = 2*(lastNll - result->
minNll());
189 log <<
"degree = " << degree
190 <<
" -log L("<<degree-1<<
") = " << lastNll
191 <<
" -log L(" << degree <<
") = " << result->
minNll()
193 <<
" P(chi^2_1 > q) = " <<
TMath::Prob(q,1) << endl;;
197 lastNll = result->
minNll();
202 log <<
"------ End Bernstein Correction Log --------" << endl;
213 const char* nominalName,
215 const char* dataName,
217 TH1F* samplingDistExtra,
225 if (!x || !nominal || !data) {
226 cout <<
"Error: wrong name for pdf or variable or dataset ! " << std::endl;
231 std::stringstream
log;
232 log <<
"------ Begin Bernstein Correction Log --------" << endl;
238 vector<RooRealVar*> coefficients;
241 for(
int i = 0; i<=degree+1; ++i) {
243 std::stringstream str;
247 "Bernstein basis poly coefficient",
251 if(i<degree) coeffNull.
add(*newCoef);
252 if(i<=degree) coeff.
add(*newCoef);
253 coeffExtra.
add(*newCoef);
254 coefficients.push_back(newCoef);
259 =
new RooBernstein(
"poly",
"Bernstein poly", *x, coeff);
263 =
new RooBernstein(
"polyNull",
"Bernstein poly", *x, coeffNull);
267 =
new RooBernstein(
"polyExtra",
"Bernstein poly", *x, coeffExtra);
271 =
new RooEffProd(
"corrected",
"",*nominal,*poly);
274 =
new RooEffProd(
"correctedNull",
"",*nominal,*polyNull);
277 =
new RooEffProd(
"correctedExtra",
"",*nominal,*polyExtra);
280 cout <<
"made pdfs, make toy generator" << endl;
290 if (printLevel < 0) {
298 for(
int i=0; i<nToys; ++i){
299 cout <<
"on toy " << i << endl;
323 samplingDist->
Fill(q);
324 samplingDistExtra->
Fill(qExtra);
326 cout <<
"NLL Results: null " << resultNull->
minNll() <<
" ref = " << result->
minNll() <<
" extra" << resultExtra->
minNll() << endl;
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Minos(Bool_t flag=kTRUE)
Bernstein basis polynomials are positive-definite in the range [0,1].
RooFit::MsgLevel globalKillBelow() const
tomato 1-D histogram with a float per channel (see TH1 documentation)}
static RooMsgService & instance()
Return reference to singleton instance.
Double_t Prob(Double_t chi2, Int_t ndf)
Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf...
RooDataSet is a container class to hold N-dimensional binned data.
Int_t ImportCorrectedPdf(RooWorkspace *, const char *, const char *, const char *)
Main method for Bernstein correction.
RooRealVar represents a fundamental (non-derived) real valued object.
RooHistPdf implements a probablity density function sampled from a multidimensional histogram...
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found...
virtual void setVal(Double_t value)
Set value of variable to 'value'.
static const std::string & DefaultMinimizerType()
void setConstant(Bool_t value=kTRUE)
RooCmdArg Minimizer(const char *type, const char *alg=0)
static int DefaultPrintLevel()
void setGlobalKillBelow(RooFit::MsgLevel level)
RooAbsData is the common abstract base class for binned and unbinned datasets.
RooDataSet is a container class to hold unbinned data.
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.
Namespace for the RooStats classes.
RooCmdArg Hesse(Bool_t flag=kTRUE)
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
RooCmdArg Save(Bool_t flag=kTRUE)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
BernsteinCorrection is a utility in RooStats to augment a nominal PDF with a polynomial correction te...
The RooWorkspace is a persistable container for RooFit projects.
virtual Int_t numEntries() const