Logo ROOT   6.16/01
Reference Guide
rf103_interprfuncs.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'BASIC FUNCTIONALITY' RooFit tutorial macro #103
5///
6/// Interpreted functions and p.d.f.s
7///
8/// \macro_image
9/// \macro_output
10/// \macro_code
11/// \author 07/2008 - Wouter Verkerke
12
13
14#include "RooRealVar.h"
15#include "RooDataSet.h"
16#include "RooGaussian.h"
17#include "TCanvas.h"
18#include "TAxis.h"
19#include "RooPlot.h"
20#include "RooFitResult.h"
21#include "RooGenericPdf.h"
22#include "RooConstVar.h"
23
24using namespace RooFit ;
25
26
28{
29 // ----------------------------------------------------
30 // G e n e r i c i n t e r p r e t e d p . d . f .
31 // ====================================================
32
33 // Declare observable x
34 RooRealVar x("x","x",-20,20) ;
35
36 // C o n s t r u c t g e n e r i c p d f f r o m i n t e r p r e t e d e x p r e s s i o n
37 // -------------------------------------------------------------------------------------------------
38
39 // To construct a proper p.d.f, the formula expression is explicitly normalized internally by dividing
40 // it by a numeric integral of the expression over x in the range [-20,20]
41 //
42 RooRealVar alpha("alpha","alpha",5,0.1,10) ;
43 RooGenericPdf genpdf("genpdf","genpdf","(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))",RooArgSet(x,alpha)) ;
44
45
46 // S a m p l e , f i t a n d p l o t g e n e r i c p d f
47 // ---------------------------------------------------------------
48
49 // Generate a toy dataset from the interpreted p.d.f
50 RooDataSet* data = genpdf.generate(x,10000) ;
51
52 // Fit the interpreted p.d.f to the generated data
53 genpdf.fitTo(*data) ;
54
55 // Make a plot of the data and the p.d.f overlaid
56 RooPlot* xframe = x.frame(Title("Interpreted expression pdf")) ;
57 data->plotOn(xframe) ;
58 genpdf.plotOn(xframe) ;
59
60
61 // -----------------------------------------------------------------------------------------------------------
62 // S t a n d a r d p . d . f a d j u s t w i t h i n t e r p r e t e d h e l p e r f u n c t i o n
63 // ==========================================================================================================
64 // Make a gauss(x,sqrt(mean2),sigma) from a standard RooGaussian
65
66
67 // C o n s t r u c t s t a n d a r d p d f w i t h f o r m u l a r e p l a c i n g p a r a m e t e r
68 // ------------------------------------------------------------------------------------------------------------
69
70 // Construct parameter mean2 and sigma
71 RooRealVar mean2("mean2","mean^2",10,0,200) ;
72 RooRealVar sigma("sigma","sigma",3,0.1,10) ;
73
74 // Construct interpreted function mean = sqrt(mean^2)
75 RooFormulaVar mean("mean","mean","sqrt(mean2)",mean2) ;
76
77 // Construct a gaussian g2(x,sqrt(mean2),sigma) ;
78 RooGaussian g2("g2","h2",x,mean,sigma) ;
79
80
81 // G e n e r a t e t o y d a t a
82 // ---------------------------------
83
84 // Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian dataset with mean 10 and width 3
85 RooGaussian g1("g1","g1",x,RooConst(10),RooConst(3)) ;
86 RooDataSet* data2 = g1.generate(x,1000) ;
87
88
89 // F i t a n d p l o t t a i l o r e d s t a n d a r d p d f
90 // -------------------------------------------------------------------
91
92 // Fit g2 to data from g1
93 RooFitResult* r = g2.fitTo(*data2,Save()) ;
94 r->Print() ;
95
96 // Plot data on frame and overlay projection of g2
97 RooPlot* xframe2 = x.frame(Title("Tailored Gaussian pdf")) ;
98 data2->plotOn(xframe2) ;
99 g2.plotOn(xframe2) ;
100
101
102 // Draw all frames on a canvas
103 TCanvas* c = new TCanvas("rf103_interprfuncs","rf103_interprfuncs",800,400) ;
104 c->Divide(2) ;
105 c->cd(1) ; gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;
106 c->cd(2) ; gPad->SetLeftMargin(0.15) ; xframe2->GetYaxis()->SetTitleOffset(1.4) ; xframe2->Draw() ;
107
108
109}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define c(i)
Definition: RSha256.hxx:101
#define gPad
Definition: TVirtualPad.h:286
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
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:531
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:550
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
RooConstVar & RooConst(Double_t val)
RooCmdArg Save(Bool_t flag=kTRUE)
const char * Title
Definition: TXMLSetup.cxx:67