Logo ROOT  
Reference Guide
rf313_paramranges.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 ///
5 /// Multidimensional models: working with parametrized ranges to define non-rectangular
6 /// regions for fitting and integration
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 ///
12 /// \date 07/2008
13 /// \author Wouter Verkerke
14 
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "RooConstVar.h"
19 #include "RooPolynomial.h"
20 #include "RooProdPdf.h"
21 #include "TCanvas.h"
22 #include "TAxis.h"
23 #include "RooPlot.h"
24 using namespace RooFit;
25 
26 void rf313_paramranges()
27 {
28 
29  // C r e a t e 3 D p d f
30  // -------------------------
31 
32  // Define observable (x,y,z)
33  RooRealVar x("x", "x", 0, 10);
34  RooRealVar y("y", "y", 0, 10);
35  RooRealVar z("z", "z", 0, 10);
36 
37  // Define 3 dimensional pdf
38  RooRealVar z0("z0", "z0", -0.1, 1);
39  RooPolynomial px("px", "px", x, RooConst(0));
40  RooPolynomial py("py", "py", y, RooConst(0));
41  RooPolynomial pz("pz", "pz", z, z0);
42  RooProdPdf pxyz("pxyz", "pxyz", RooArgSet(px, py, pz));
43 
44  // D e f i n e d n o n - r e c t a n g u l a r r e g i o n R i n ( x , y , z )
45  // -------------------------------------------------------------------------------------
46 
47  //
48  // R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
49  //
50 
51  // Construct range parametrized in "R" in y [ 0.1*x, 0.9*x ]
52  RooFormulaVar ylo("ylo", "0.1*x", x);
53  RooFormulaVar yhi("yhi", "0.9*x", x);
54  y.setRange("R", ylo, yhi);
55 
56  // Construct parametrized ranged "R" in z [ 0, 0.1*y^2 ]
57  RooFormulaVar zlo("zlo", "0.0*y", y);
58  RooFormulaVar zhi("zhi", "0.1*y*y", y);
59  z.setRange("R", zlo, zhi);
60 
61  // C a l c u l a t e i n t e g r a l o f n o r m a l i z e d p d f i n R
62  // ----------------------------------------------------------------------------------
63 
64  // Create integral over normalized pdf model over x,y,z in "R" region
65  RooAbsReal *intPdf = pxyz.createIntegral(RooArgSet(x, y, z), RooArgSet(x, y, z), "R");
66 
67  // Plot value of integral as function of pdf parameter z0
68  RooPlot *frame = z0.frame(Title("Integral of pxyz over x,y,z in region R"));
69  intPdf->plotOn(frame);
70 
71  new TCanvas("rf313_paramranges", "rf313_paramranges", 600, 600);
72  gPad->SetLeftMargin(0.15);
73  frame->GetYaxis()->SetTitleOffset(1.6);
74  frame->Draw();
75 
76  return;
77 }
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
RooAbsReal::createIntegral
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, 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
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:561
RooGaussian.h
x
Double_t x[n]
Definition: legend1.C:17
RooAbsReal
Definition: RooAbsReal.h:61
TCanvas.h
RooDataSet.h
RooPolynomial.h
RooPlot::frame
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition: RooPlot.cxx:249
RooFormulaVar
Definition: RooFormulaVar.h:30
RooProdPdf.h
RooFit
Definition: RooCFunction1Binding.h:29
RooPolynomial
Definition: RooPolynomial.h:28
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1256
RooPlot
Definition: RooPlot.h:44
rf313_paramranges
Definition: rf313_paramranges.py:1
y
Double_t y[n]
Definition: legend1.C:17
RooRealVar.h
RooConstVar.h
RooAbsReal::plotOn
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
Definition: RooAbsReal.cxx:1714
TCanvas
Definition: TCanvas.h:23
TAxis.h
gPad
#define gPad
Definition: TVirtualPad.h:287
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
RooFit::Title
RooCmdArg Title(const char *name)
Definition: RooGlobalFunc.cxx:173
RooArgSet
Definition: RooArgSet.h:28
RooFit::RooConst
RooConstVar & RooConst(Double_t val)
Definition: RooGlobalFunc.cxx:341