Logo ROOT  
Reference Guide
rf110_normintegration.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Basic functionality: normalization and integration of pdfs, construction of cumulative distribution
5 /// monodimensional functions
6 ///
7 /// \macro_image
8 /// \macro_output
9 /// \macro_code
10 ///
11 /// \date July 2008
12 /// \author Wouter Verkerke
13 
14 #include "RooRealVar.h"
15 #include "RooGaussian.h"
16 #include "RooConstVar.h"
17 #include "RooAbsReal.h"
18 #include "RooPlot.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 using namespace RooFit;
22 
24 {
25  // S e t u p m o d e l
26  // ---------------------
27 
28  // Create observables x,y
29  RooRealVar x("x", "x", -10, 10);
30 
31  // Create pdf gaussx(x,-2,3)
32  RooGaussian gx("gx", "gx", x, RooConst(-2), RooConst(3));
33 
34  // R e t r i e v e r a w & n o r m a l i z e d v a l u e s o f R o o F i t p . d . f . s
35  // --------------------------------------------------------------------------------------------------
36 
37  // Return 'raw' unnormalized value of gx
38  cout << "gx = " << gx.getVal() << endl;
39 
40  // Return value of gx normalized over x in range [-10,10]
41  RooArgSet nset(x);
42  cout << "gx_Norm[x] = " << gx.getVal(&nset) << endl;
43 
44  // Create object representing integral over gx
45  // which is used to calculate gx_Norm[x] == gx / gx_Int[x]
46  RooAbsReal *igx = gx.createIntegral(x);
47  cout << "gx_Int[x] = " << igx->getVal() << endl;
48 
49  // I n t e g r a t e n o r m a l i z e d p d f o v e r s u b r a n g e
50  // ----------------------------------------------------------------------------
51 
52  // Define a range named "signal" in x from -5,5
53  x.setRange("signal", -5, 5);
54 
55  // Create an integral of gx_Norm[x] over x in range "signal"
56  // This is the fraction of of pdf gx_Norm[x] which is in the
57  // range named "signal"
58  RooAbsReal *igx_sig = gx.createIntegral(x, NormSet(x), Range("signal"));
59  cout << "gx_Int[x|signal]_Norm[x] = " << igx_sig->getVal() << endl;
60 
61  // C o n s t r u c t c u m u l a t i v e d i s t r i b u t i o n f u n c t i o n f r o m p d f
62  // -----------------------------------------------------------------------------------------------------
63 
64  // Create the cumulative distribution function of gx
65  // i.e. calculate Int[-10,x] gx(x') dx'
66  RooAbsReal *gx_cdf = gx.createCdf(x);
67 
68  // Plot cdf of gx versus x
69  RooPlot *frame = x.frame(Title("cdf of Gaussian pdf"));
70  gx_cdf->plotOn(frame);
71 
72  // Draw plot on canvas
73  new TCanvas("rf110_normintegration", "rf110_normintegration", 600, 600);
74  gPad->SetLeftMargin(0.15);
75  frame->GetYaxis()->SetTitleOffset(1.6);
76  frame->Draw();
77 }
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
RooAbsReal.h
RooAbsReal::createIntegral
RooAbsReal * createIntegral(const RooArgSet &iset, 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
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:574
RooGaussian.h
RooAbsReal::getVal
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:91
x
Double_t x[n]
Definition: legend1.C:17
RooGaussian
Plain Gaussian p.d.f.
Definition: RooGaussian.h:24
RooAbsReal
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:61
TCanvas.h
RooFit
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition: RooCFunction1Binding.h:29
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1258
RooPlot
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
RooRealVar.h
rf110_normintegration
Definition: rf110_normintegration.py:1
RooConstVar.h
RooAbsReal::plotOn
virtual RooPlot * plotOn(RooPlot *frame, 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(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
Definition: RooAbsReal.cxx:1730
RooFit::NormSet
RooCmdArg NormSet(const RooArgSet &nset)
Definition: RooGlobalFunc.cxx:262
TCanvas
The Canvas class.
Definition: TCanvas.h:23
TAxis.h
RooFit::Range
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Definition: RooGlobalFunc.cxx:53
gPad
#define gPad
Definition: TVirtualPad.h:287
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:37
RooFit::Title
RooCmdArg Title(const char *name)
Definition: RooGlobalFunc.cxx:176
RooArgSet
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:29
RooFit::RooConst
RooConstVar & RooConst(Double_t val)
Definition: RooGlobalFunc.cxx:347