ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs_bernsteinCorrection.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'Bernstein Correction' RooStats tutorial macro
4 // author: Kyle Cranmer
5 // date March. 2009
6 //
7 // This tutorial shows usage of a the BernsteinCorrection utility in RooStats.
8 // The idea is that one has a distribution coming either from data or Monte Carlo
9 // (called "reality" in the macro) and a nominal model that is not sufficiently
10 // flexible to take into account the real distribution. One wants to take into
11 // account the systematic associated with this imperfect modeling by augmenting
12 // the nominal model with some correction term (in this case a polynomial).
13 // The BernsteinCorrection utility will import into your workspace a corrected model
14 // given by nominal(x) * poly_N(x), where poly_N is an n-th order polynomial in
15 // the Bernstein basis. The degree N of the polynomial is chosen by specifying the tolerance
16 // one has in adding an extra term to the polynomial.
17 // The Bernstein basis is nice because it only has positive-definite terms
18 // and works well with PDFs.
19 // Finally, the macro makes a plot of:
20 // - the data (drawn from 'reality'),
21 // - the best fit of the nominal model (blue)
22 // - and the best fit corrected model.
23 /////////////////////////////////////////////////////////////////////////
24 
25 #ifndef __CINT__
26 #include "RooGlobalFunc.h"
27 #endif
28 #include "RooDataSet.h"
29 #include "RooRealVar.h"
30 #include "RooConstVar.h"
31 #include "RooBernstein.h"
32 #include "TCanvas.h"
33 #include "RooAbsPdf.h"
34 #include "RooFit.h"
35 #include "RooFitResult.h"
36 #include "RooPlot.h"
37 #include <string>
38 #include <vector>
39 #include <stdio.h>
40 #include <sstream>
41 #include <iostream>
42 
43 #include "RooProdPdf.h"
44 #include "RooAddPdf.h"
45 #include "RooGaussian.h"
46 #include "RooNLLVar.h"
47 #include "RooMinuit.h"
48 #include "RooProfileLL.h"
49 #include "RooWorkspace.h"
50 
52 
53 // use this order for safety on library loading
54 using namespace RooFit;
55 using namespace RooStats;
56 
57 
58 //____________________________________
60 
61  // set range of observable
62  Double_t lowRange = -1, highRange =5;
63 
64  // make a RooRealVar for the observable
65  RooRealVar x("x", "x", lowRange, highRange);
66 
67  // true model
68  RooGaussian narrow("narrow","",x,RooConst(0.), RooConst(.8));
69  RooGaussian wide("wide","",x,RooConst(0.), RooConst(2.));
70  RooAddPdf reality("reality","",RooArgList(narrow, wide), RooConst(0.8));
71 
72  RooDataSet* data = reality.generate(x,1000);
73 
74  // nominal model
75  RooRealVar sigma("sigma","",1.,0,10);
76  RooGaussian nominal("nominal","",x,RooConst(0.), sigma);
77 
78  RooWorkspace* wks = new RooWorkspace("myWorksspace");
79 
80  wks->import(*data, Rename("data"));
81  wks->import(nominal);
82 
83  // use Minuit2
85 
86  // The tolerance sets the probability to add an unnecessary term.
87  // lower tolerance will add fewer terms, while higher tolerance
88  // will add more terms and provide a more flexible function.
89  Double_t tolerance = 0.05;
90  BernsteinCorrection bernsteinCorrection(tolerance);
91  Int_t degree = bernsteinCorrection.ImportCorrectedPdf(wks,"nominal","x","data");
92 
93  if (degree < 0) {
94  Error("rs_bernsteinCorrection","Bernstein correction failed ! ");
95  return;
96  }
97 
98  cout << " Correction based on Bernstein Poly of degree " << degree << endl;
99 
100 
101  RooPlot* frame = x.frame();
102  data->plotOn(frame);
103  // plot the best fit nominal model in blue
105  nominal.fitTo(*data,PrintLevel(0),Minimizer(minimType));
106  nominal.plotOn(frame);
107 
108  // plot the best fit corrected model in red
109  RooAbsPdf* corrected = wks->pdf("corrected");
110  if (!corrected) return;
111 
112  // fit corrected model
113  corrected->fitTo(*data,PrintLevel(0),Minimizer(minimType) );
114  corrected->plotOn(frame,LineColor(kRed));
115 
116  // plot the correction term (* norm constant) in dashed green
117  // should make norm constant just be 1, not depend on binning of data
118  RooAbsPdf* poly = wks->pdf("poly");
119  if (poly)
120  poly->plotOn(frame,LineColor(kGreen), LineStyle(kDashed));
121 
122  // this is a switch to check the sampling distribution
123  // of -2 log LR for two comparisons:
124  // the first is for n-1 vs. n degree polynomial corrections
125  // the second is for n vs. n+1 degree polynomial corrections
126  // Here we choose n to be the one chosen by the tolerance
127  // critereon above, eg. n = "degree" in the code.
128  // Setting this to true is takes about 10 min.
129  bool checkSamplingDist = true;
130  int numToyMC = 20; // increse this value for sensible results
131 
132  TCanvas* c1 = new TCanvas();
133  if(checkSamplingDist) {
134  c1->Divide(1,2);
135  c1->cd(1);
136  }
137  frame->Draw();
138  gPad->Update();
139 
140  if(checkSamplingDist) {
141  // check sampling dist
143  TH1F* samplingDist = new TH1F("samplingDist","",20,0,10);
144  TH1F* samplingDistExtra = new TH1F("samplingDistExtra","",20,0,10);
145  bernsteinCorrection.CreateQSamplingDist(wks,"nominal","x","data",samplingDist, samplingDistExtra, degree,numToyMC);
146 
147  c1->cd(2);
148  samplingDistExtra->SetLineColor(kRed);
149  samplingDistExtra->Draw();
150  samplingDist->Draw("same");
151  }
152 }
153 
RooCmdArg LineColor(Color_t color)
TCanvas * c1
Definition: legend1.C:2
Definition: Rtypes.h:61
RooCmdArg PrintLevel(Int_t code)
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
Definition: Rtypes.h:61
static void SetDefaultPrintLevel(int level)
Double_t x[n]
Definition: legend1.C:17
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooCmdArg LineStyle(Style_t style)
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
Int_t ImportCorrectedPdf(RooWorkspace *, const char *, const char *, const char *)
Main method for Bernstein correction.
void rs_bernsteinCorrection()
const Double_t sigma
virtual RooPlot * plotOn(RooPlot *frame, 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()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:626
void Error(const char *location, const char *msgfmt,...)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
static const std::string & DefaultMinimizerType()
virtual RooPlot * plotOn(RooPlot *frame, 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(), const RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Plot (project) PDF on specified frame.
Definition: RooAbsPdf.h:105
RooCmdArg Minimizer(const char *type, const char *alg=0)
RooCmdArg Rename(const char *suffix)
RooPlot * frame(const RooCmdArg &arg1, 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()) const
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
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.
double Double_t
Definition: RtypesCore.h:55
RooConstVar & RooConst(Double_t val)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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...
Definition: RooAbsPdf.cxx:1702
#define gPad
Definition: TVirtualPad.h:288
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.
Definition: RooAbsPdf.cxx:1056
BernsteinCorrection is a utility in RooStats to augment a nominal PDF with a polynomial correction te...
static void SetDefaultMinimizer(const char *type, const char *algo=0)
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559