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 p.d.fs, construction of cumulative distribution
5/// monodimensional functions
6///
7/// \macro_image
8/// \macro_output
9/// \macro_code
10/// \author 07/2008 - Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooGaussian.h"
14#include "RooConstVar.h"
15#include "RooAbsReal.h"
16#include "RooPlot.h"
17#include "TCanvas.h"
18#include "TAxis.h"
19using namespace RooFit;
20
22{
23 // S e t u p m o d e l
24 // ---------------------
25
26 // Create observables x,y
27 RooRealVar x("x", "x", -10, 10);
28
29 // Create p.d.f. gaussx(x,-2,3)
30 RooGaussian gx("gx", "gx", x, RooConst(-2), RooConst(3));
31
32 // 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
33 // --------------------------------------------------------------------------------------------------
34
35 // Return 'raw' unnormalized value of gx
36 cout << "gx = " << gx.getVal() << endl;
37
38 // Return value of gx normalized over x in range [-10,10]
39 RooArgSet nset(x);
40 cout << "gx_Norm[x] = " << gx.getVal(&nset) << endl;
41
42 // Create object representing integral over gx
43 // which is used to calculate gx_Norm[x] == gx / gx_Int[x]
44 RooAbsReal *igx = gx.createIntegral(x);
45 cout << "gx_Int[x] = " << igx->getVal() << endl;
46
47 // 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
48 // ----------------------------------------------------------------------------
49
50 // Define a range named "signal" in x from -5,5
51 x.setRange("signal", -5, 5);
52
53 // Create an integral of gx_Norm[x] over x in range "signal"
54 // This is the fraction of of p.d.f. gx_Norm[x] which is in the
55 // range named "signal"
56 RooAbsReal *igx_sig = gx.createIntegral(x, NormSet(x), Range("signal"));
57 cout << "gx_Int[x|signal]_Norm[x] = " << igx_sig->getVal() << endl;
58
59 // 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
60 // -----------------------------------------------------------------------------------------------------
61
62 // Create the cumulative distribution function of gx
63 // i.e. calculate Int[-10,x] gx(x') dx'
64 RooAbsReal *gx_cdf = gx.createCdf(x);
65
66 // Plot cdf of gx versus x
67 RooPlot *frame = x.frame(Title("c.d.f of Gaussian p.d.f"));
68 gx_cdf->plotOn(frame);
69
70 // Draw plot on canvas
71 new TCanvas("rf110_normintegration", "rf110_normintegration", 600, 600);
72 gPad->SetLeftMargin(0.15);
73 frame->GetYaxis()->SetTitleOffset(1.6);
74 frame->Draw();
75}
#define gPad
Definition: TVirtualPad.h:286
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
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.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
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:562
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg NormSet(const RooArgSet &nset)
RooConstVar & RooConst(Double_t val)
const char * Title
Definition: TXMLSetup.cxx:67
Ta Range(0, 0, 1, 1)