Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 pdf with a convolution operator in terms of a parameter
5///
6/// This tutorial requires FFT3 to be enabled.
7///
8/// \macro_image
9/// \macro_code
10/// \macro_output
11///
12/// \date April 2009
13/// \author Wouter Verkerke
14
15#include "RooRealVar.h"
16#include "RooDataHist.h"
17#include "RooGaussian.h"
18#include "RooGenericPdf.h"
19#include "RooFormulaVar.h"
20#include "RooFFTConvPdf.h"
21#include "RooPlot.h"
22#include "TCanvas.h"
23#include "TAxis.h"
24#include "TH2.h"
25using namespace RooFit;
26
27void rf211_paramconv()
28{
29 // S e t u p c o m p o n e n t p d f s
30 // ---------------------------------------
31
32 // Gaussian g(x ; mean,sigma)
33 RooRealVar x("x", "x", -10, 10);
34 RooRealVar mean("mean", "mean", -3, 3);
35 RooRealVar sigma("sigma", "sigma", 0.5, 0.1, 10);
36 RooGaussian modelx("gx", "gx", x, mean, sigma);
37
38 // Block function in mean
39 RooRealVar a("a", "a", 2, 1, 10);
40 RooGenericPdf model_mean("model_mean", "abs(mean)<a", RooArgList(mean, a));
41
42 // Convolution in mean parameter model = g(x,mean,sigma) (x) block(mean)
43 x.setBins(1000, "cache");
44 mean.setBins(50, "cache");
45 RooFFTConvPdf model("model", "model", mean, modelx, model_mean);
46
47 // Configure convolution to construct a 2-D cache in (x,mean)
48 // rather than a 1-d cache in mean that needs to be recalculated
49 // for each value of x
50 model.setCacheObservables(x);
51 model.setBufferFraction(1.0);
52
53 // Integrate model over mean projModel = Int model dmean
54 RooAbsPdf *projModel = model.createProjection(mean);
55
56 // Generate 1000 toy events
57 std::unique_ptr<RooDataHist> d{projModel->generateBinned(x, 1000)};
58
59 // Fit pdf to toy data
60 projModel->fitTo(*d, Verbose(), PrintLevel(-1));
61
62 // Plot data and fitted pdf
63 RooPlot *frame = x.frame(Bins(25));
64 d->plotOn(frame);
65 projModel->plotOn(frame);
66
67 // Make 2d histogram of model(x;mean)
68 TH1 *hh = model.createHistogram("hh", x, Binning(50), YVar(mean, Binning(50)), ConditionalObservables(mean));
69 hh->SetTitle("histogram of model(x|mean)");
71
72 // Draw frame on canvas
73 TCanvas *c = new TCanvas("rf211_paramconv", "rf211_paramconv", 800, 400);
74 c->Divide(2);
75 c->cd(1);
76 gPad->SetLeftMargin(0.15);
77 frame->GetYaxis()->SetTitleOffset(1.4);
78 frame->Draw();
79 c->cd(2);
80 gPad->SetLeftMargin(0.20);
81 hh->GetZaxis()->SetTitleOffset(2.5);
82 hh->Draw("surf");
83}
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
@ kBlue
Definition Rtypes.h:66
#define gPad
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:124
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:157
virtual RooFit::OwningPtr< RooDataHist > generateBinned(const RooArgSet &whatVars, double nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}) const
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition RooAbsPdf.h:110
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
PDF for the numerical (FFT) convolution of two PDFs.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Implementation of a probability density function that takes a RooArgList of servers and a C++ express...
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
TAxis * GetZaxis()
Definition TH1.h:326
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6686
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg Bins(Int_t nbin)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Verbose(bool flag=true)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
RooCmdArg Binning(const RooAbsBinning &binning)
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...
Definition JSONIO.h:26