Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf211_paramconv.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Addition and convolution: working with a pdf with a convolution operator in terms of a parameter

This tutorial requires FFT3 to be enabled.

#include "RooRealVar.h"
#include "RooDataHist.h"
#include "RooGaussian.h"
#include "RooGenericPdf.h"
#include "RooFormulaVar.h"
#include "RooFFTConvPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH2.h"
using namespace RooFit;
{
// S e t u p c o m p o n e n t p d f s
// ---------------------------------------
// Gaussian g(x ; mean,sigma)
RooRealVar x("x", "x", -10, 10);
RooRealVar mean("mean", "mean", -3, 3);
RooRealVar sigma("sigma", "sigma", 0.5, 0.1, 10);
RooGaussian modelx("gx", "gx", x, mean, sigma);
// Block function in mean
RooRealVar a("a", "a", 2, 1, 10);
RooGenericPdf model_mean("model_mean", "abs(mean)<a", RooArgList(mean, a));
// Convolution in mean parameter model = g(x,mean,sigma) (x) block(mean)
x.setBins(1000, "cache");
mean.setBins(50, "cache");
RooFFTConvPdf model("model", "model", mean, modelx, model_mean);
// Configure convolution to construct a 2-D cache in (x,mean)
// rather than a 1-d cache in mean that needs to be recalculated
// for each value of x
model.setCacheObservables(x);
model.setBufferFraction(1.0);
// Integrate model over mean projModel = Int model dmean
RooAbsPdf *projModel = model.createProjection(mean);
// Generate 1000 toy events
RooDataHist *d = projModel->generateBinned(x, 1000);
// Fit pdf to toy data
projModel->fitTo(*d, Verbose(), PrintLevel(-1));
// Plot data and fitted pdf
RooPlot *frame = x.frame(Bins(25));
d->plotOn(frame);
projModel->plotOn(frame);
// Make 2d histogram of model(x;mean)
TH1 *hh = model.createHistogram("hh", x, Binning(50), YVar(mean, Binning(50)), ConditionalObservables(mean));
hh->SetTitle("histogram of model(x|mean)");
// Draw frame on canvas
TCanvas *c = new TCanvas("rf211_paramconv", "rf211_paramconv", 800, 400);
c->Divide(2);
c->cd(1);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.4);
frame->Draw();
c->cd(2);
gPad->SetLeftMargin(0.20);
hh->GetZaxis()->SetTitleOffset(2.5);
hh->Draw("surf");
}
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
PrintLevel
Definition RooMinuit.h:6
@ kBlue
Definition Rtypes.h:66
#define gPad
virtual RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, const RooLinkedList &cmdList={})
Fit PDF to given dataset.
virtual RooFit::OwningPtr< RooDataHist > generateBinned(const RooArgSet &whatVars, double 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:112
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
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 override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:126
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
PDF for the numerical (FFT) convolution of two PDFs.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
A RooPlot is a 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:239
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
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
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:58
TAxis * GetZaxis()
Definition TH1.h:324
void SetTitle(const char *title) override
Change (i.e.
Definition TH1.cxx:6700
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
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 Common.h:18
[#1] INFO:Eval -- RooRealVar::setRange(mean) new range named 'refrange_fft_model' created with bounds [-3,3]
[#0] WARNING:Eval -- The FFT convolution 'model' will run with 50 bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number of bins of the observable 'mean'.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x564b05461000 with pdf gx_CONV_model_mean_CACHE_Obs[mean,x]_NORM_mean for nset (mean) with code 0
[#0] WARNING:Eval -- The FFT convolution 'model' will run with 50 bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number of bins of the observable 'mean'.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x564b056569b0 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean]_NORM_x_mean for nset (x,mean) with code 1
[#0] WARNING:Eval -- The FFT convolution 'model' will run with 50 bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number of bins of the observable 'mean'.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x564b05664a40 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean]_NORM_x_mean for nset (x,mean) with code 1 from preexisting content.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for a: using 0.5
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for sigma: using 0.2
prevFCN = 2171.275755 a=2.012, sigma=0.5, [#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2.012,
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2, sigma=0.5047,
prevFCN = 2171.286528 sigma=0.4953,
prevFCN = 2171.267762 sigma=0.5029,
prevFCN = 2171.281998 sigma=0.4971,
prevFCN = 2171.270547 a=2.012, sigma=0.5,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2.012,
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2, sigma=0.5029,
prevFCN = 2171.281998 sigma=0.4971,
prevFCN = 2171.270547 a=2.013, sigma=0.4843,
prevFCN = 2171.259881 a=2.009, sigma=0.4884,
prevFCN = 2171.260992 a=2.025, sigma=0.4843,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.134,
prevFCN = 2171.692149 a=1.898,
prevFCN = 2172.378568 a=2.025,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.013, sigma=0.4871,
prevFCN = 2171.26042 sigma=0.4815,
prevFCN = 2171.260367 sigma=0.4843,
prevFCN = 2171.259881 a=2.025,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.134,
prevFCN = 2171.692149 a=1.898,
prevFCN = 2172.378568 a=2.013, sigma=0.4871,
prevFCN = 2171.26042 sigma=0.4815,
prevFCN = 2171.260367 a=2.015, sigma=0.4843,
prevFCN = 2171.259881 a=2.011,
prevFCN = 2171.259881 a=2.013, sigma=0.4848,
prevFCN = 2171.259907 sigma=0.4837,
prevFCN = 2171.259896 a=2.134, sigma=0.4871,
prevFCN = 2171.720718 a=2.065, sigma=0.4512,
prevFCN = 2171.332894 a=2.03, sigma=0.473,
prevFCN = 2171.26812 a=2.02, sigma=0.4794,
prevFCN = 2171.261381 a=2.016, sigma=0.482,
prevFCN = 2171.260198 a=2.014, sigma=0.4832,
prevFCN = 2171.25995 a=2.014, sigma=0.4837,
prevFCN = 2171.259895 a=2.013, sigma=0.484,
prevFCN = 2171.259883 a=2.013, sigma=0.4841,
prevFCN = 2171.259881 a=2.013, sigma=0.4842,
prevFCN = 2171.25988 a=2.025,
prevFCN = 2171.25988 a=2.001,
prevFCN = 2171.25988 a=2.134,
prevFCN = 2171.691427 a=1.898,
prevFCN = 2172.379556 a=2.025,
prevFCN = 2171.25988 a=2.001,
prevFCN = 2171.25988 a=2.013, sigma=0.487,
prevFCN = 2171.260398 sigma=0.4814,
prevFCN = 2171.260398 sigma=0.4842,
prevFCN = 2171.25988 a=2.025,
prevFCN = 2171.25988 a=2.001,
prevFCN = 2171.25988 a=2.134,
prevFCN = 2171.691427 a=1.898,
prevFCN = 2172.379556 a=2.013, sigma=0.487,
prevFCN = 2171.260398 sigma=0.4814,
prevFCN = 2171.260398 a=2.015, sigma=0.4842,
prevFCN = 2171.25988 a=2.011,
prevFCN = 2171.25988 a=2.013, sigma=0.4848,
prevFCN = 2171.259901 sigma=0.4836,
prevFCN = 2171.259901 sigma=0.4843,
prevFCN = 2171.259881 sigma=0.4841,
prevFCN = 2171.259881 a=2.134, sigma=0.487,
prevFCN = 2171.720107 a=2.013, sigma=0.4842, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#0] WARNING:Eval -- The FFT convolution 'model' will run with 50 bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number of bins of the observable 'mean'.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x564b05666010 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean]_NORM_x_mean for nset (x,mean) with code 1
[#0] WARNING:Eval -- The FFT convolution 'model' will run with 50 bins. A decent accuracy for difficult convolutions is typically only reached with n >= 1000. Suggest to increase the number of bins of the observable 'mean'.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x564b05666010 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean]_NORM_x for nset (x) with code 1 from preexisting content.
Date
April 2009
Author
Wouter Verkerke

Definition in file rf211_paramconv.C.