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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: making 2/3 dimensional plots of pdfs and datasets

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooProdPdf.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
#include "RooPlot.h"
using namespace RooFit;
{
// C r e a t e 2 D m o d e l a n d d a t a s e t
// -----------------------------------------------------
// Create observables
RooRealVar x("x", "x", -5, 5);
RooRealVar y("y", "y", -5, 5);
// Create parameters
RooRealVar a0("a0", "a0", -3.5, -5, 5);
RooRealVar a1("a1", "a1", -1.5, -1, 1);
RooRealVar sigma("sigma", "width of gaussian", 1.5);
// Create interpreted function f(y) = a0 - a1*sqrt(10*abs(y))
RooFormulaVar fy("fy", "a0-a1*sqrt(10*abs(y))", RooArgSet(y, a0, a1));
// Create gauss(x,f(y),s)
RooGaussian model("model", "Gaussian with shifting mean", x, fy, sigma);
// Sample dataset from gauss(x,y)
std::unique_ptr<RooDataSet> data{model.generate({x, y}, 10000)};
// M a k e 2 D p l o t s o f d a t a a n d m o d e l
// -------------------------------------------------------------
// Create and fill ROOT 2D histogram (20x20 bins) with contents of dataset
TH1 *hh_data = data->createHistogram("x,y", Binning(20), Binning(20));
// Create and fill ROOT 2D histogram (50x50 bins) with sampling of pdf
// TH2D* hh_pdf = model.createHistogram("hh_model",x,Binning(50),YVar(y,Binning(50))) ;
TH1 *hh_pdf = model.createHistogram("x,y", 50, 50);
hh_pdf->SetLineColor(kBlue);
// C r e a t e 3 D m o d e l a n d d a t a s e t
// -----------------------------------------------------
// Create observables
RooRealVar z("z", "z", -5, 5);
RooGaussian gz("gz", "gz", z, 0.0, 2.0);
RooProdPdf model3("model3", "model3", RooArgSet(model, gz));
std::unique_ptr<RooDataSet> data3{model3.generate({x, y, z}, 10000)};
// M a k e 3 D p l o t s o f d a t a a n d m o d e l
// -------------------------------------------------------------
// Create and fill ROOT 2D histogram (8x8x8 bins) with contents of dataset
TH1 *hh_data3 = data3->createHistogram("hh_data3", x, Binning(8), YVar(y, Binning(8)), ZVar(z, Binning(8)));
// Create and fill ROOT 2D histogram (20x20x20 bins) with sampling of pdf
TH1 *hh_pdf3 = model3.createHistogram("hh_model3", x, Binning(20), YVar(y, Binning(20)), ZVar(z, Binning(20)));
hh_pdf3->SetFillColor(kBlue);
TCanvas *c1 = new TCanvas("rf309_2dimplot", "rf309_2dimplot", 800, 800);
c1->Divide(2, 2);
c1->cd(1);
gPad->SetLeftMargin(0.15);
hh_data->GetZaxis()->SetTitleOffset(1.4);
hh_data->Draw("lego");
c1->cd(2);
gPad->SetLeftMargin(0.20);
hh_pdf->GetZaxis()->SetTitleOffset(2.5);
hh_pdf->Draw("surf");
c1->cd(3);
gPad->SetLeftMargin(0.15);
hh_data->GetZaxis()->SetTitleOffset(1.4);
hh_data->Draw("box");
c1->cd(4);
gPad->SetLeftMargin(0.15);
hh_pdf->GetZaxis()->SetTitleOffset(2.5);
hh_pdf->Draw("cont3");
TCanvas *c2 = new TCanvas("rf309_3dimplot", "rf309_3dimplot", 800, 400);
c2->Divide(2);
c2->cd(1);
gPad->SetLeftMargin(0.15);
hh_data3->GetZaxis()->SetTitleOffset(1.4);
hh_data3->Draw("lego");
c2->cd(2);
gPad->SetLeftMargin(0.15);
hh_pdf3->GetZaxis()->SetTitleOffset(1.4);
hh_pdf3->Draw("iso");
}
@ 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
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:39
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 SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
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:327
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
RooCmdArg ZVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg Binning(const RooAbsBinning &binning)
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
[#0] WARNING:InputArguments -- The parameter 'sigma' with range [-inf, inf] of the RooGaussian 'model' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooIntegrator1D to calculate Int(y)
Date
July 2008
Author
Wouter Verkerke

Definition in file rf309_ndimplot.C.