Logo ROOT   6.10/09
Reference Guide
rf110_normintegration.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #110
4 //
5 // Examples on normalization of p.d.f.s,
6 // integration of p.d.fs, construction
7 // of cumulative distribution functions from p.d.f.s
8 // in one dimension
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooGaussian.h"
19 #include "RooAbsReal.h"
20 #include "RooPlot.h"
21 #include "TCanvas.h"
22 using namespace RooFit ;
23 
24 class TestBasic110 : public RooFitTestUnit
25 {
26 public:
27  TestBasic110(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Normalization of p.d.f.s in 1D",refFile,writeRef,verbose) {} ;
28  Bool_t testCode() {
29 
30  // S e t u p m o d e l
31  // ---------------------
32 
33  // Create observables x,y
34  RooRealVar x("x","x",-10,10) ;
35 
36  // Create p.d.f. gaussx(x,-2,3)
37  RooGaussian gx("gx","gx",x,RooConst(-2),RooConst(3)) ;
38 
39 
40  // 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
41  // --------------------------------------------------------------------------------------------------
42 
43  // Return 'raw' unnormalized value of gx
44  regValue(gx.getVal(),"rf110_gx") ;
45 
46  // Return value of gx normalized over x in range [-10,10]
47  RooArgSet nset(x) ;
48 
49  regValue(gx.getVal(&nset),"rf110_gx_Norm[x]") ;
50 
51  // Create object representing integral over gx
52  // which is used to calculate gx_Norm[x] == gx / gx_Int[x]
53  RooAbsReal* igx = gx.createIntegral(x) ;
54  regValue(igx->getVal(),"rf110_gx_Int[x]") ;
55 
56 
57  // 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
58  // ----------------------------------------------------------------------------
59 
60  // Define a range named "signal" in x from -5,5
61  x.setRange("signal",-5,5) ;
62 
63  // Create an integral of gx_Norm[x] over x in range "signal"
64  // This is the fraction of of p.d.f. gx_Norm[x] which is in the
65  // range named "signal"
66  RooAbsReal* igx_sig = gx.createIntegral(x,NormSet(x),Range("signal")) ;
67  regValue(igx_sig->getVal(),"rf110_gx_Int[x|signal]_Norm[x]") ;
68 
69 
70 
71  // 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
72  // -----------------------------------------------------------------------------------------------------
73 
74  // Create the cumulative distribution function of gx
75  // i.e. calculate Int[-10,x] gx(x') dx'
76  RooAbsReal* gx_cdf = gx.createCdf(x) ;
77 
78  // Plot cdf of gx versus x
79  RooPlot* frame = x.frame(Title("c.d.f of Gaussian p.d.f")) ;
80  gx_cdf->plotOn(frame) ;
81 
82 
83  regPlot(frame,"rf110_plot1") ;
84 
85  delete igx ;
86  delete igx_sig ;
87  delete gx_cdf ;
88 
89  return kTRUE ;
90  }
91 } ;
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
RooCmdArg NormSet(const RooArgSet &nset)
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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:501
bool verbose
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooConstVar & RooConst(Double_t val)
const Bool_t kTRUE
Definition: RtypesCore.h:91