Logo ROOT   6.16/01
Reference Guide
rf211_paramconv.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #211
5///
6/// Working a with a p.d.f. with a convolution operator in terms
7/// of a parameter
8///
9/// (require ROOT to be compiled with --enable-fftw3)
10///
11/// \macro_image
12/// \macro_output
13/// \macro_code
14/// \author 04/2009 - Wouter Verkerke
15
16
17#include "RooRealVar.h"
18#include "RooDataHist.h"
19#include "RooGaussian.h"
20#include "RooGenericPdf.h"
21#include "RooFormulaVar.h"
22#include "RooFFTConvPdf.h"
23#include "RooPlot.h"
24#include "TCanvas.h"
25#include "TAxis.h"
26#include "TH2.h"
27using namespace RooFit ;
28
29
30void rf211_paramconv()
31{
32 // S e t u p c o m p o n e n t p d f s
33 // ---------------------------------------
34
35 // Gaussian g(x ; mean,sigma)
36 RooRealVar x("x","x",-10,10) ;
37 RooRealVar mean("mean","mean",-3,3) ;
38 RooRealVar sigma("sigma","sigma",0.5,0.1,10) ;
39 RooGaussian modelx("gx","gx",x,mean,sigma) ;
40
41 // Block function in mean
42 RooRealVar a("a","a",2,1,10) ;
43 RooGenericPdf model_mean("model_mean","abs(mean)<a",RooArgList(mean,a)) ;
44
45 // Convolution in mean parameter model = g(x,mean,sigma) (x) block(mean)
46 x.setBins(1000,"cache") ;
47 mean.setBins(50,"cache") ;
48 RooFFTConvPdf model("model","model",mean,modelx,model_mean) ;
49
50 // Configure convolution to construct a 2-D cache in (x,mean)
51 // rather than a 1-d cache in mean that needs to be recalculated
52 // for each value of x
53 model.setCacheObservables(x) ;
54 model.setBufferFraction(1.0) ;
55
56 // Integrate model over mean projModel = Int model dmean
57 RooAbsPdf* projModel = model.createProjection(mean) ;
58
59 // Generate 1000 toy events
60 RooDataHist* d = projModel->generateBinned(x,1000) ;
61
62 // Fit p.d.f. to toy data
63 projModel->fitTo(*d,Verbose()) ;
64
65 // Plot data and fitted p.d.f.
66 RooPlot* frame = x.frame(Bins(25)) ;
67 d->plotOn(frame) ;
68 projModel->plotOn(frame) ;
69
70 // Make 2d histogram of model(x;mean)
71 TH1* hh = model.createHistogram("hh",x,Binning(50),YVar(mean,Binning(50)),ConditionalObservables(mean)) ;
72 hh->SetTitle("histogram of model(x|mean)") ;
73 hh->SetLineColor(kBlue) ;
74
75 // Draw frame on canvas
76 TCanvas* c = new TCanvas("rf211_paramconv","rf211_paramconv",800,400) ;
77 c->Divide(2) ;
78 c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
79 c->cd(2) ; gPad->SetLeftMargin(0.20) ; hh->GetZaxis()->SetTitleOffset(2.5) ; hh->Draw("surf") ;
80
81}
#define d(i)
Definition: RSha256.hxx:102
#define c(i)
Definition: RSha256.hxx:101
@ kBlue
Definition: Rtypes.h:63
#define gPad
Definition: TVirtualPad.h:286
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual RooFitResult * fitTo(RooAbsData &data, 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())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1081
virtual RooDataHist * generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition: RooAbsPdf.h:105
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 RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition: RooAbsPdf.h:119
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
PDF for the numerical (FFT) convolution of two PDFs.
Definition: RooFFTConvPdf.h:26
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
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:31
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6217
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg ConditionalObservables(const RooArgSet &set)
RooCmdArg Verbose(Bool_t flag=kTRUE)
RooCmdArg Bins(Int_t nbin)
auto * a
Definition: textangle.C:12