ROOT
v6-30
Reference Guide
Loading...
Searching...
No Matches
rf103_interprfuncs.C File Reference
Tutorials
»
RooFit Tutorials
Detailed Description
Basic functionality: interpreted functions and PDFs.
#include "
RooRealVar.h
"
#include "
RooDataSet.h
"
#include "
RooGaussian.h
"
#include "
TCanvas.h
"
#include "
TAxis.h
"
#include "
RooPlot.h
"
#include "
RooFitResult.h
"
#include "
RooGenericPdf.h
"
using namespace
RooFit
;
void
rf103_interprfuncs
()
{
// ----------------------------------------------------
// G e n e r i c i n t e r p r e t e d p . d . f .
// ====================================================
// Declare observable x
RooRealVar
x
(
"x"
,
"x"
, -20, 20);
// C o n s t r u c t g e n e r i c p d f f r o m i n t e r p r e t e d e x p r e s s i o n
// -------------------------------------------------------------------------------------------------
// To construct a proper pdf, the formula expression is explicitly normalized internally by dividing
// it by a numeric integral of the expression over x in the range [-20,20]
//
RooRealVar
alpha(
"alpha"
,
"alpha"
, 5, 0.1, 10);
RooGenericPdf
genpdf(
"genpdf"
,
"genpdf"
,
"(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))"
,
RooArgSet
(
x
, alpha));
// S a m p l e , f i t a n d p l o t g e n e r i c p d f
// ---------------------------------------------------------------
// Generate a toy dataset from the interpreted pdf
std::unique_ptr<RooDataSet>
data
{genpdf.generate(
x
, 10000)};
// Fit the interpreted pdf to the generated data
genpdf.fitTo(*
data
, PrintLevel(-1));
// Make a plot of the data and the pdf overlaid
RooPlot
*xframe =
x
.
frame
(Title(
"Interpreted expression pdf"
));
data
->plotOn(xframe);
genpdf.plotOn(xframe);
// -----------------------------------------------------------------------------------------------------------
// S t a n d a r d p . d . f a d j u s t w i t h i n t e r p r e t e d h e l p e r f u n c t i o n
// ==========================================================================================================
// Make a gauss(x,sqrt(mean2),sigma) from a standard RooGaussian
// C o n s t r u c t s t a n d a r d p d f w i t h f o r m u l a r e p l a c i n g p a r a m e t e r
// ------------------------------------------------------------------------------------------------------------
// Construct parameter mean2 and sigma
RooRealVar
mean2(
"mean2"
,
"mean^2"
, 10, 0, 200);
RooRealVar
sigma
(
"sigma"
,
"sigma"
, 3, 0.1, 10);
// Construct interpreted function mean = sqrt(mean^2)
RooFormulaVar
mean(
"mean"
,
"mean"
,
"sqrt(mean2)"
, mean2);
// Construct a gaussian g2(x,sqrt(mean2),sigma) ;
RooGaussian
g2(
"g2"
,
"h2"
,
x
, mean,
sigma
);
// G e n e r a t e t o y d a t a
// ---------------------------------
// Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian dataset with mean 10 and width 3
RooGaussian
g1(
"g1"
,
"g1"
,
x
, 10.0, 3.0);
std::unique_ptr<RooDataSet> data2{g1.generate(
x
, 1000)};
// F i t a n d p l o t t a i l o r e d s t a n d a r d p d f
// -------------------------------------------------------------------
// Fit g2 to data from g1
std::unique_ptr<RooFitResult> fitResult{g2.fitTo(*data2,
Save
(),
PrintLevel
(-1))};
fitResult->Print();
// Plot data on frame and overlay projection of g2
RooPlot
*xframe2 =
x
.
frame
(
Title
(
"Tailored Gaussian pdf"
));
data2->plotOn(xframe2);
g2.plotOn(xframe2);
// Draw all frames on a canvas
TCanvas
*
c
=
new
TCanvas
(
"rf103_interprfuncs"
,
"rf103_interprfuncs"
, 800, 400);
c
->Divide(2);
c
->cd(1);
gPad
->SetLeftMargin(0.15);
xframe->
GetYaxis
()->
SetTitleOffset
(1.4);
xframe->
Draw
();
c
->cd(2);
gPad
->SetLeftMargin(0.15);
xframe2->
GetYaxis
()->
SetTitleOffset
(1.4);
xframe2->
Draw
();
}
c
#define c(i)
Definition
RSha256.hxx:101
RooDataSet.h
RooFitResult.h
RooGaussian.h
RooGenericPdf.h
RooPlot.h
RooRealVar.h
TAxis.h
TCanvas.h
data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Definition
TGWin32VirtualXProxy.cxx:104
gPad
#define gPad
Definition
TVirtualPad.h:305
RooArgSet
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition
RooArgSet.h:55
RooFormulaVar
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition
RooFormulaVar.h:30
RooGaussian
Plain Gaussian p.d.f.
Definition
RooGaussian.h:24
RooGenericPdf
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition
RooGenericPdf.h:25
RooPlot
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition
RooPlot.h:43
RooPlot::frame
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
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition
RooPlot.cxx:1279
RooPlot::Draw
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition
RooPlot.cxx:652
RooRealVar
RooRealVar represents a variable that can be changed from the outside.
Definition
RooRealVar.h:37
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition
TAttAxis.cxx:298
TCanvas
The Canvas class.
Definition
TCanvas.h:23
RooFit::Save
RooCmdArg Save(bool flag=true)
Definition
RooGlobalFunc.cxx:564
RooFit::PrintLevel
RooCmdArg PrintLevel(Int_t code)
Definition
RooGlobalFunc.cxx:572
sigma
const Double_t sigma
Definition
h1analysisProxy.h:11
x
Double_t x[n]
Definition
legend1.C:17
RooFit
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition
JSONIO.h:26
rf103_interprfuncs
Definition
rf103_interprfuncs.py:1
xmlio::Title
const char * Title
Definition
TXMLSetup.cxx:68
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooRombergIntegrator to calculate Int(x)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooRombergIntegrator to calculate Int(x)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooRombergIntegrator to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooRombergIntegrator to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 2551.39, estimated distance to minimum: 4.39288e-06
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
mean2 1.0010e+02 +/- 1.98e+00
sigma 3.1172e+00 +/- 7.12e-02
Date
July 2008
Author
Wouter Verkerke
Definition in file
rf103_interprfuncs.C
.
tutorials
roofit
rf103_interprfuncs.C
ROOT v6-30 - Reference Guide Generated on Sun Dec 1 2024 05:18:44 (GVA Time) using Doxygen 1.9.8