ROOT
6.30.04
Reference Guide
All
Namespaces
Files
Pages
rf313_paramranges.C
Go to the documentation of this file.
1
/// \file
2
/// \ingroup tutorial_roofit
3
/// \notebook -js
4
/// Multidimensional models: working with parametrized ranges to define non-rectangular regions for fitting and
5
/// integration
6
///
7
/// \macro_image
8
/// \macro_output
9
/// \macro_code
10
/// \author 07/2008 - Wouter Verkerke
11
12
#include "
RooRealVar.h
"
13
#include "
RooDataSet.h
"
14
#include "
RooGaussian.h
"
15
#include "
RooConstVar.h
"
16
#include "
RooPolynomial.h
"
17
#include "
RooProdPdf.h
"
18
#include "
TCanvas.h
"
19
#include "
TAxis.h
"
20
#include "
RooPlot.h
"
21
using namespace
RooFit;
22
23
void
rf313_paramranges()
24
{
25
26
// C r e a t e 3 D p d f
27
// -------------------------
28
29
// Define observable (x,y,z)
30
RooRealVar x(
"x"
,
"x"
, 0, 10);
31
RooRealVar y(
"y"
,
"y"
, 0, 10);
32
RooRealVar z(
"z"
,
"z"
, 0, 10);
33
34
// Define 3 dimensional pdf
35
RooRealVar z0(
"z0"
,
"z0"
, -0.1, 1);
36
RooPolynomial px(
"px"
,
"px"
, x, RooConst(0));
37
RooPolynomial py(
"py"
,
"py"
, y, RooConst(0));
38
RooPolynomial pz(
"pz"
,
"pz"
, z, z0);
39
RooProdPdf pxyz(
"pxyz"
,
"pxyz"
, RooArgSet(px, py, pz));
40
41
// 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 )
42
// -------------------------------------------------------------------------------------
43
44
//
45
// R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
46
//
47
48
// Construct range parametrized in "R" in y [ 0.1*x, 0.9*x ]
49
RooFormulaVar ylo(
"ylo"
,
"0.1*x"
, x);
50
RooFormulaVar yhi(
"yhi"
,
"0.9*x"
, x);
51
y.setRange(
"R"
, ylo, yhi);
52
53
// Construct parametrized ranged "R" in z [ 0, 0.1*y^2 ]
54
RooFormulaVar zlo(
"zlo"
,
"0.0*y"
, y);
55
RooFormulaVar zhi(
"zhi"
,
"0.1*y*y"
, y);
56
z.setRange(
"R"
, zlo, zhi);
57
58
// 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
59
// ----------------------------------------------------------------------------------
60
61
// Create integral over normalized pdf model over x,y,z in "R" region
62
RooAbsReal *intPdf = pxyz.createIntegral(RooArgSet(x, y, z), RooArgSet(x, y, z),
"R"
);
63
64
// Plot value of integral as function of pdf parameter z0
65
RooPlot *frame = z0.frame(Title(
"Integral of pxyz over x,y,z in region R"
));
66
intPdf->plotOn(frame);
67
68
new
TCanvas(
"rf313_paramranges"
,
"rf313_paramranges"
, 600, 600);
69
gPad->SetLeftMargin(0.15);
70
frame->GetYaxis()->SetTitleOffset(1.6);
71
frame->Draw();
72
73
return
;
74
}
RooProdPdf.h
RooRealVar.h
RooGaussian.h
TCanvas.h
RooConstVar.h
RooPlot.h
RooPolynomial.h
TAxis.h
RooDataSet.h
tutorials
roofit
rf313_paramranges.C
Generated on Tue May 5 2020 14:03:49 for ROOT by
1.8.5