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

Detailed Description

View in nbviewer Open in SWAN
Organization and simultaneous fits: working with named parameter sets and parameter snapshots in workspaces

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooChebychev.h"
#include "RooAddPdf.h"
#include "RooWorkspace.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TFile.h"
#include "TH1.h"
using namespace RooFit;
void fillWorkspace(RooWorkspace &w);
{
// C r e a t e m o d e l a n d d a t a s e t
// -----------------------------------------------
fillWorkspace(*w);
// Exploit convention encoded in named set "parameters" and "observables"
// to use workspace contents w/o need for introspected
RooAbsPdf *model = w->pdf("model");
// Generate data from pdf in given observables
std::unique_ptr<RooDataSet> data{model->generate(*w->set("observables"), 1000)};
// Fit model to data
model->fitTo(*data, PrintLevel(-1));
// Plot fitted model and data on frame of first (only) observable
RooPlot *frame = ((RooRealVar *)w->set("observables")->first())->frame();
data->plotOn(frame);
model->plotOn(frame);
// Overlay plot with model with reference parameters as stored in snapshots
w->loadSnapshot("reference_fit");
model->plotOn(frame, LineColor(kRed));
w->loadSnapshot("reference_fit_bkgonly");
model->plotOn(frame, LineColor(kRed), LineStyle(kDashed));
// Draw the frame on the canvas
new TCanvas("rf510_wsnamedsets", "rf503_wsnamedsets", 600, 600);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.4);
frame->Draw();
// Print workspace contents
w->Print();
// Workspace will remain in memory after macro finishes
gDirectory->Add(w);
}
void fillWorkspace(RooWorkspace &w)
{
// C r e a t e m o d e l
// -----------------------
// Declare observable x
RooRealVar x("x", "x", 0, 10);
// Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
RooRealVar mean("mean", "mean of gaussians", 5, 0, 10);
RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
RooRealVar sigma2("sigma2", "width of gaussians", 1);
RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
// Build Chebychev polynomial pdf
RooRealVar a0("a0", "a0", 0.5, 0., 1.);
RooRealVar a1("a1", "a1", 0.2, 0., 1.);
RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
// Sum the signal components into a composite signal pdf
RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
// Sum the composite signal and background
RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
// Import model into pdf
w.import(model);
// E n c o d e d e f i n i t i o n o f p a r a m e t e r s i n w o r k s p a c e
// ---------------------------------------------------------------------------------------
// Define named sets "parameters" and "observables", which list which variables should be considered
// parameters and observables by the users convention
//
// Variables appearing in sets _must_ live in the workspace already, or the autoImport flag
// of defineSet must be set to import them on the fly. Named sets contain only references
// to the original variables, therefore the value of observables in named sets already
// reflect their 'current' value
std::unique_ptr<RooArgSet> params{model.getParameters(x)};
w.defineSet("parameters", *params);
w.defineSet("observables", x);
// E n c o d e r e f e r e n c e v a l u e f o r p a r a m e t e r s i n w o r k s p a c e
// ---------------------------------------------------------------------------------------------------
// Define a parameter 'snapshot' in the pdf
// Unlike a named set, a parameter snapshot stores an independent set of values for
// a given set of variables in the workspace. The values can be stored and reloaded
// into the workspace variable objects using the loadSnapshot() and saveSnapshot()
// methods. A snapshot saves the value of each variable, any errors that are stored
// with it as well as the 'Constant' flag that is used in fits to determine if a
// parameter is kept fixed or not.
// Do a dummy fit to a (supposedly) reference dataset here and store the results
// of that fit into a snapshot
std::unique_ptr<RooDataSet> refData{model.generate(x, 10000)};
model.fitTo(*refData, PrintLevel(-1));
// The true flag imports the values of the objects in (*params) into the workspace
// If not set, the present values of the workspace parameters objects are stored
w.saveSnapshot("reference_fit", *params, true);
// Make another fit with the signal component forced to zero
// and save those parameters too
bkgfrac.setVal(1);
bkgfrac.setConstant(true);
bkgfrac.removeError();
model.fitTo(*refData, PrintLevel(-1));
w.saveSnapshot("reference_fit_bkgonly", *params, true);
}
@ kRed
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
#define gDirectory
Definition TDirectory.h:384
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:123
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:156
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:57
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Chebychev polynomial p.d.f.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
TAxis * GetYaxis() const
Definition RooPlot.cxx:1279
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:652
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
Persistable container for RooFit projects.
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg PrintLevel(Int_t code)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-inf, inf] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-inf, inf] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooChebychev::bkg
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::x
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a0
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::bkgfrac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::sig
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::mean
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sig1frac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig2
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma2
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooWorkspace(w) w contents
variables
---------
(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2,x)
p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 1
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 1/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 0.999388/1
RooGaussian::sig1[ x=x mean=mean sigma=sigma1 ] = 0.999291
RooGaussian::sig2[ x=x mean=mean sigma=sigma2 ] = 0.999823
parameter snapshots
-------------------
reference_fit = (a0=0.500613 +/- 0.023199,a1=0.160315 +/- 0.0373121,bkgfrac=0.504699 +/- 0.0113933,mean=5.01883 +/- 0.0101222,sigma1=0.5[C],sig1frac=0.8179 +/- 0.0374037,sigma2=1[C])
reference_fit_bkgonly = (a0=0.474264 +/- 0,a1=6.8252e-12 +/- 0,bkgfrac=1[C],mean=5.01883 +/- 0,sigma1=0.5[C],sig1frac=0.8179 +/- 0,sigma2=1[C])
named sets
----------
observables:(x)
parameters:(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2)
Date
April 2009
Author
Wouter Verkerke

Definition in file rf510_wsnamedsets.C.