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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: marginizalization of multi-dimensional pdfs through integration

#include "RooRealVar.h"
#include "RooDataHist.h"
#include "RooGaussian.h"
#include "RooProdPdf.h"
#include "RooPolyVar.h"
#include "TH1.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;
{
// C r e a t e p d f m ( x , y ) = g x ( x | y ) * g ( y )
// --------------------------------------------------------------
// Increase default precision of numeric integration
// as this exercise has high sensitivity to numeric integration precision
// Create observables
RooRealVar x("x", "x", -5, 5);
RooRealVar y("y", "y", -2, 2);
// Create function f(y) = a0 + a1*y
RooRealVar a0("a0", "a0", 0);
RooRealVar a1("a1", "a1", -1.5, -3, 1);
RooPolyVar fy("fy", "fy", y, RooArgSet(a0, a1));
// Create gaussx(x,f(y),sx)
RooRealVar sigmax("sigmax", "width of gaussian", 0.5);
RooGaussian gaussx("gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax);
// Create gaussy(y,0,2)
RooGaussian gaussy("gaussy", "Gaussian in y", y, 0.0, 2.0);
// Create gaussx(x,sx|y) * gaussy(y)
RooProdPdf model("model", "gaussx(x|y)*gaussy(y)", gaussy, Conditional(gaussx, x));
// M a r g i n a l i z e m ( x , y ) t o m ( x )
// ----------------------------------------------------
// modelx(x) = Int model(x,y) dy
RooAbsPdf *modelx = model.createProjection(y);
// U s e m a r g i n a l i z e d p . d . f . a s r e g u l a r 1 - D p . d . f .
// ------------------------------------------------------------------------------------------
// Sample 1000 events from modelx
std::unique_ptr<RooAbsData> data{modelx->generateBinned(x, 1000)};
// Fit modelx to toy data
modelx->fitTo(*data, Verbose(), PrintLevel(-1));
// Plot modelx over data
RooPlot *frame = x.frame(40);
data->plotOn(frame);
modelx->plotOn(frame);
// Make 2D histogram of model(x,y)
TH1 *hh = model.createHistogram("x,y");
TCanvas *c = new TCanvas("rf315_projectpdf", "rf315_projectpdf", 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 c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#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:123
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:156
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:109
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
void setEpsRel(double newEpsRel)
Set relative convergence criteria (convergence if std::abs(Err)/abs(Int)<newEpsRel)
void setEpsAbs(double newEpsAbs)
Set absolute convergence criteria (convergence if std::abs(Err)<newEpsAbs)
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:1279
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:652
Class RooPolyVar is a RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficie...
Definition RooPolyVar.h:25
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
RooRealVar represents a 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:58
TAxis * GetZaxis()
Definition TH1.h:324
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3067
Double_t y[n]
Definition legend1.C:17
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
[#0] WARNING:InputArguments -- The parameter 'sigmax' with range [-inf, inf] of the RooGaussian 'gaussx' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (gaussx)
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for a1: using 0.4
prevFCN = 1900.156536 a1=-1.488,
prevFCN = 1899.96972 a1=-1.512,
prevFCN = 1900.566064 a1=-1.499,
prevFCN = 1900.127678 a1=-1.501,
prevFCN = 1900.187622 a1=-1.484,
prevFCN = 1899.958651 a1=-1.483,
prevFCN = 1899.959674 a1=-1.484,
prevFCN = 1899.958586 a1=-1.484,
prevFCN = 1899.958651 a1=-1.483,
prevFCN = 1899.959674 a1=-1.484,
prevFCN = 1899.958586 a1=-1.483,
prevFCN = 1899.958779 a1=-1.484,
prevFCN = 1899.958562 a1=-1.484,
prevFCN = 1899.958651 a1=-1.483,
prevFCN = 1899.958779 a1=-1.484,
prevFCN = 1899.958562 a1=-1.484,
prevFCN = 1899.958674 a1=-1.484,
prevFCN = 1899.95863 a1=-1.484, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
Date
July 2008
Author
Wouter Verkerke

Definition in file rf315_projectpdf.C.