Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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_code
9/// \macro_output
10///
11/// \date July 2008
12/// \author Wouter Verkerke
13
14#include "RooRealVar.h"
15#include "RooGaussian.h"
16#include "RooAbsReal.h"
17#include "RooPlot.h"
18#include "TCanvas.h"
19#include "TAxis.h"
20using namespace RooFit;
21
23{
24 // S e t u p m o d e l
25 // ---------------------
26
27 // Create observables x,y
28 RooRealVar x("x", "x", -10, 10);
29
30 // Create pdf gaussx(x,-2,3)
31 RooGaussian gx("gx", "gx", x, -2.0, 3.0);
32
33 // 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
34 // --------------------------------------------------------------------------------------------------
35
36 // Return 'raw' unnormalized value of gx
37 cout << "gx = " << gx.getVal() << endl;
38
39 // Return value of gx normalized over x in range [-10,10]
40 RooArgSet nset(x);
41 cout << "gx_Norm[x] = " << gx.getVal(&nset) << endl;
42
43 // Create object representing integral over gx
44 // which is used to calculate gx_Norm[x] == gx / gx_Int[x]
45 std::unique_ptr<RooAbsReal> igx{gx.createIntegral(x)};
46 cout << "gx_Int[x] = " << igx->getVal() << endl;
47
48 // 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
49 // ----------------------------------------------------------------------------
50
51 // Define a range named "signal" in x from -5,5
52 x.setRange("signal", -5, 5);
53
54 // Create an integral of gx_Norm[x] over x in range "signal"
55 // This is the fraction of of pdf gx_Norm[x] which is in the
56 // range named "signal"
57 std::unique_ptr<RooAbsReal> igx_sig{gx.createIntegral(x, NormSet(x), Range("signal"))};
58 cout << "gx_Int[x|signal]_Norm[x] = " << igx_sig->getVal() << endl;
59
60 // 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
61 // -----------------------------------------------------------------------------------------------------
62
63 // Create the cumulative distribution function of gx
64 // i.e. calculate Int[-10,x] gx(x') dx'
65 std::unique_ptr<RooAbsReal> gx_cdf{gx.createCdf(x)};
66
67 // Plot cdf of gx versus x
68 RooPlot *frame = x.frame(Title("cdf of Gaussian pdf"));
69 gx_cdf->plotOn(frame);
70
71 // Draw plot on canvas
72 new TCanvas("rf110_normintegration", "rf110_normintegration", 600, 600);
73 gPad->SetLeftMargin(0.15);
74 frame->GetYaxis()->SetTitleOffset(1.6);
75 frame->Draw();
76}
#define gPad
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg NormSet(Args_t &&... argsOrArgSet)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68
Ta Range(0, 0, 1, 1)