Logo ROOT   6.10/09
Reference Guide
rf103_interprfuncs.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #103
4 //
5 // Interpreted functions and p.d.f.s
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "TCanvas.h"
20 #include "RooPlot.h"
21 #include "RooFitResult.h"
22 #include "RooGenericPdf.h"
23 
24 using namespace RooFit ;
25 
26 
27 class TestBasic103 : public RooFitTestUnit
28 {
29 public:
30  TestBasic103(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Interpreted expression p.d.f.",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  /////////////////////////////////////////////////////////
34  // G e n e r i c i n t e r p r e t e d p . d . f . //
35  /////////////////////////////////////////////////////////
36 
37  // Declare observable x
38  RooRealVar x("x","x",-20,20) ;
39 
40  // 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
41  // -------------------------------------------------------------------------------------------------
42 
43  // To construct a proper p.d.f, the formula expression is explicitly normalized internally by dividing
44  // it by a numeric integral of the expresssion over x in the range [-20,20]
45  //
46  RooRealVar alpha("alpha","alpha",5,0.1,10) ;
47  RooGenericPdf genpdf("genpdf","genpdf","(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))",RooArgSet(x,alpha)) ;
48 
49 
50  // 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
51  // ---------------------------------------------------------------
52 
53  // Generate a toy dataset from the interpreted p.d.f
54  RooDataSet* data = genpdf.generate(x,10000) ;
55 
56  // Fit the interpreted p.d.f to the generated data
57  genpdf.fitTo(*data) ;
58 
59  // Make a plot of the data and the p.d.f overlaid
60  RooPlot* xframe = x.frame(Title("Interpreted expression pdf")) ;
61  data->plotOn(xframe) ;
62  genpdf.plotOn(xframe) ;
63 
64 
65  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
66  // 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 //
67  // //
68  // Make a gauss(x,sqrt(mean2),sigma) from a standard RooGaussian //
69  // //
70  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////
71 
72 
73  // 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
74  // ------------------------------------------------------------------------------------------------------------
75 
76  // Construct parameter mean2 and sigma
77  RooRealVar mean2("mean2","mean^2",10,0,200) ;
78  RooRealVar sigma("sigma","sigma",3,0.1,10) ;
79 
80  // Construct interpreted function mean = sqrt(mean^2)
81  RooFormulaVar mean("mean","mean","sqrt(mean2)",mean2) ;
82 
83  // Construct a gaussian g2(x,sqrt(mean2),sigma) ;
84  RooGaussian g2("g2","h2",x,mean,sigma) ;
85 
86 
87  // G e n e r a t e t o y d a t a
88  // ---------------------------------
89 
90  // Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian dataset with mean 10 and width 3
91  RooGaussian g1("g1","g1",x,RooConst(10),RooConst(3)) ;
92  RooDataSet* data2 = g1.generate(x,1000) ;
93 
94 
95  // 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
96  // -------------------------------------------------------------------
97 
98  // Fit g2 to data from g1
99  RooFitResult* r = g2.fitTo(*data2,Save()) ;
100 
101  // Plot data on frame and overlay projection of g2
102  RooPlot* xframe2 = x.frame(Title("Tailored Gaussian pdf")) ;
103  data2->plotOn(xframe2) ;
104  g2.plotOn(xframe2) ;
105 
106  regPlot(xframe,"rf103_plot1") ;
107  regPlot(xframe2,"rf103_plot2") ;
108  regResult(r,"rf103_fit1") ;
109 
110  delete data ;
111  delete data2 ;
112 
113  return kTRUE ;
114  }
115 } ;
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
const Double_t sigma
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
TRandom2 r(17)
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg Save(Bool_t flag=kTRUE)
RooConstVar & RooConst(Double_t val)
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
const Bool_t kTRUE
Definition: RtypesCore.h:91