Logo ROOT  
Reference Guide
rf211_paramconv.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// Addition and convolution: working with a p.d.f. with a convolution operator in terms of a parameter
5///
6/// This tutorial requires FFT3 to be enabled.
7///
8/// \macro_image
9/// \macro_output
10/// \macro_code
11/// \author 04/2009 - Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataHist.h"
15#include "RooGaussian.h"
16#include "RooGenericPdf.h"
17#include "RooFormulaVar.h"
18#include "RooFFTConvPdf.h"
19#include "RooPlot.h"
20#include "TCanvas.h"
21#include "TAxis.h"
22#include "TH2.h"
23using namespace RooFit;
24
25void rf211_paramconv()
26{
27 // S e t u p c o m p o n e n t p d f s
28 // ---------------------------------------
29
30 // Gaussian g(x ; mean,sigma)
31 RooRealVar x("x", "x", -10, 10);
32 RooRealVar mean("mean", "mean", -3, 3);
33 RooRealVar sigma("sigma", "sigma", 0.5, 0.1, 10);
34 RooGaussian modelx("gx", "gx", x, mean, sigma);
35
36 // Block function in mean
37 RooRealVar a("a", "a", 2, 1, 10);
38 RooGenericPdf model_mean("model_mean", "abs(mean)<a", RooArgList(mean, a));
39
40 // Convolution in mean parameter model = g(x,mean,sigma) (x) block(mean)
41 x.setBins(1000, "cache");
42 mean.setBins(50, "cache");
43 RooFFTConvPdf model("model", "model", mean, modelx, model_mean);
44
45 // Configure convolution to construct a 2-D cache in (x,mean)
46 // rather than a 1-d cache in mean that needs to be recalculated
47 // for each value of x
48 model.setCacheObservables(x);
49 model.setBufferFraction(1.0);
50
51 // Integrate model over mean projModel = Int model dmean
52 RooAbsPdf *projModel = model.createProjection(mean);
53
54 // Generate 1000 toy events
55 RooDataHist *d = projModel->generateBinned(x, 1000);
56
57 // Fit p.d.f. to toy data
58 projModel->fitTo(*d, Verbose());
59
60 // Plot data and fitted p.d.f.
61 RooPlot *frame = x.frame(Bins(25));
62 d->plotOn(frame);
63 projModel->plotOn(frame);
64
65 // Make 2d histogram of model(x;mean)
66 TH1 *hh = model.createHistogram("hh", x, Binning(50), YVar(mean, Binning(50)), ConditionalObservables(mean));
67 hh->SetTitle("histogram of model(x|mean)");
69
70 // Draw frame on canvas
71 TCanvas *c = new TCanvas("rf211_paramconv", "rf211_paramconv", 800, 400);
72 c->Divide(2);
73 c->cd(1);
74 gPad->SetLeftMargin(0.15);
75 frame->GetYaxis()->SetTitleOffset(1.4);
76 frame->Draw();
77 c->cd(2);
78 gPad->SetLeftMargin(0.20);
79 hh->GetZaxis()->SetTitleOffset(2.5);
80 hh->Draw("surf");
81}
#define d(i)
Definition: RSha256.hxx:102
#define c(i)
Definition: RSha256.hxx:101
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:286
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()) const
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition: RooAbsPdf.h:104
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:1286
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:3362
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:118
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
The RooDataHist 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:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
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:6333
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
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