Logo ROOT  
Reference Guide
rf315_projectpdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 ///
5 /// Multidimensional models: marginizalization of multi-dimensional p.d.f.s through integration
6 ///
7 /// \macro_image
8 /// \macro_output
9 /// \macro_code
10 ///
11 /// \date 07/2008
12 /// \author Wouter Verkerke
13 
14 #include "RooRealVar.h"
15 #include "RooDataHist.h"
16 #include "RooGaussian.h"
17 #include "RooProdPdf.h"
18 #include "RooPolyVar.h"
19 #include "TH1.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 #include "RooNumIntConfig.h"
24 #include "RooConstVar.h"
25 using namespace RooFit;
26 
27 void rf315_projectpdf()
28 {
29  // C r e a t e p d f m ( x , y ) = g x ( x | y ) * g ( y )
30  // --------------------------------------------------------------
31 
32  // Increase default precision of numeric integration
33  // as this exercise has high sensitivity to numeric integration precision
36 
37  // Create observables
38  RooRealVar x("x", "x", -5, 5);
39  RooRealVar y("y", "y", -2, 2);
40 
41  // Create function f(y) = a0 + a1*y
42  RooRealVar a0("a0", "a0", 0);
43  RooRealVar a1("a1", "a1", -1.5, -3, 1);
44  RooPolyVar fy("fy", "fy", y, RooArgSet(a0, a1));
45 
46  // Create gaussx(x,f(y),sx)
47  RooRealVar sigmax("sigmax", "width of gaussian", 0.5);
48  RooGaussian gaussx("gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax);
49 
50  // Create gaussy(y,0,2)
51  RooGaussian gaussy("gaussy", "Gaussian in y", y, RooConst(0), RooConst(2));
52 
53  // Create gaussx(x,sx|y) * gaussy(y)
54  RooProdPdf model("model", "gaussx(x|y)*gaussy(y)", gaussy, Conditional(gaussx, x));
55 
56  // M a r g i n a l i z e m ( x , y ) t o m ( x )
57  // ----------------------------------------------------
58 
59  // modelx(x) = Int model(x,y) dy
60  RooAbsPdf *modelx = model.createProjection(y);
61 
62  // 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 .
63  // ------------------------------------------------------------------------------------------
64 
65  // Sample 1000 events from modelx
66  RooAbsData *data = modelx->generateBinned(x, 1000);
67 
68  // Fit modelx to toy data
69  modelx->fitTo(*data, Verbose());
70 
71  // Plot modelx over data
72  RooPlot *frame = x.frame(40);
73  data->plotOn(frame);
74  modelx->plotOn(frame);
75 
76  // Make 2D histogram of model(x,y)
77  TH1 *hh = model.createHistogram("x,y");
78  hh->SetLineColor(kBlue);
79 
80  TCanvas *c = new TCanvas("rf315_projectpdf", "rf315_projectpdf", 800, 400);
81  c->Divide(2);
82  c->cd(1);
83  gPad->SetLeftMargin(0.15);
84  frame->GetYaxis()->SetTitleOffset(1.4);
85  frame->Draw();
86  c->cd(2);
87  gPad->SetLeftMargin(0.20);
88  hh->GetZaxis()->SetTitleOffset(2.5);
89  hh->Draw("surf");
90 }
c
#define c(i)
Definition: RSha256.hxx:119
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
e
#define e(i)
Definition: RSha256.hxx:121
RooAbsData
Definition: RooAbsData.h:46
RooNumIntConfig.h
RooGaussian.h
x
Double_t x[n]
Definition: legend1.C:17
RooGaussian
Definition: RooGaussian.h:25
rf315_projectpdf
Definition: rf315_projectpdf.py:1
TCanvas.h
TH1::GetZaxis
TAxis * GetZaxis()
Definition: TH1.h:319
RooAbsPdf::plotOn
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
RooPolyVar.h
RooProdPdf.h
RooFit
Definition: RooCFunction1Binding.h:29
RooAbsData::plotOn
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
Definition: RooAbsData.cxx:547
RooFit::Verbose
RooCmdArg Verbose(Bool_t flag=kTRUE)
Definition: RooGlobalFunc.cxx:186
RooDataHist.h
RooPlot.h
RooNumIntConfig::setEpsRel
void setEpsRel(Double_t newEpsRel)
Set relative convergence criteria (convergence if abs(Err)/abs(Int)<newEpsRel)
Definition: RooNumIntConfig.cxx:268
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1256
RooPlot
Definition: RooPlot.h:44
y
Double_t y[n]
Definition: legend1.C:17
RooRealVar.h
RooPolyVar
Definition: RooPolyVar.h:28
RooConstVar.h
RooFit::Conditional
RooCmdArg Conditional(const RooArgSet &pdfSet, const RooArgSet &depSet, Bool_t depsAreCond=kFALSE)
Definition: RooGlobalFunc.cxx:225
RooAbsReal::defaultIntegratorConfig
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Definition: RooAbsReal.cxx:3531
TCanvas
Definition: TCanvas.h:23
RooNumIntConfig::setEpsAbs
void setEpsAbs(Double_t newEpsAbs)
Set absolute convergence criteria (convergence if abs(Err)<newEpsAbs)
Definition: RooNumIntConfig.cxx:238
TAxis.h
TH1
Definition: TH1.h:57
kBlue
@ kBlue
Definition: Rtypes.h:66
gPad
#define gPad
Definition: TVirtualPad.h:287
RooAbsPdf::generateBinned
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
make_cnn_model.model
model
Definition: make_cnn_model.py:6
RooAbsPdf
Definition: RooAbsPdf.h:40
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
Definition: RooRealVar.h:35
RooProdPdf
Definition: RooProdPdf.h:33
TH1.h
RooArgSet
Definition: RooArgSet.h:28
TH1::Draw
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2997
RooFit::RooConst
RooConstVar & RooConst(Double_t val)
Definition: RooGlobalFunc.cxx:341
RooAbsPdf::fitTo
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:1261