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

␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-1e+30, 1e+30] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-1e+30, 1e+30] 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
[#0] WARNING:Eval -- Evaluating RooAddPdf without a defined normalization set. This can lead to ambiguos coefficients definition and incorrect results. Use RooAddPdf::fixCoefNormalization(nset) to provide a normalization set for defining uniquely RooAddPdf coefficients!
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#0] WARNING:Eval -- Evaluating RooAddPdf without a defined normalization set. This can lead to ambiguos coefficients definition and incorrect results. Use RooAddPdf::fixCoefNormalization(nset) to provide a normalization set for defining uniquely RooAddPdf coefficients!
[#1] INFO:Minization -- createNLL: caching constraint set under name CONSTR_OF_PDF_model_FOR_OBS_x with 0 entries
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 a0 5.00000e-01 1.00000e-01 0.00000e+00 1.00000e+00
2 a1 2.00000e-01 1.00000e-01 0.00000e+00 1.00000e+00
3 bkgfrac 5.00000e-01 1.00000e-01 0.00000e+00 1.00000e+00
4 mean 5.00000e+00 1.00000e+00 0.00000e+00 1.00000e+01
5 sig1frac 8.00000e-01 1.00000e-01 0.00000e+00 1.00000e+00
**********
** 3 **SET ERR 0.5
**********
**********
** 4 **SET PRINT 1
**********
**********
** 5 **SET STR 1
**********
NOW USING STRATEGY 1: TRY TO BALANCE SPEED AGAINST RELIABILITY
**********
** 6 **MIGRAD 2500 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=1953.02 FROM MIGRAD STATUS=INITIATE 14 CALLS 15 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a0 5.00000e-01 1.00000e-01 2.01358e-01 -5.24449e+00
2 a1 2.00000e-01 1.00000e-01 2.57889e-01 -8.31158e+00
3 bkgfrac 5.00000e-01 1.00000e-01 2.01358e-01 -1.67086e+01
4 mean 5.00000e+00 1.00000e+00 2.01358e-01 -1.20616e+02
5 sig1frac 8.00000e-01 1.00000e-01 2.57889e-01 -8.99606e+00
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=1949.67 FROM MIGRAD STATUS=CONVERGED 164 CALLS 165 TOTAL
EDM=2.46934e-05 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a0 5.35942e-01 6.82140e-02 4.16868e-03 9.78465e-04
2 a1 1.95503e-01 7.79188e-02 5.51158e-03 2.08158e-03
3 bkgfrac 5.40925e-01 2.13219e-02 1.19176e-03 1.14661e-02
4 mean 5.01914e+00 2.93156e-02 1.77988e-04 4.20703e-01
5 sig1frac 9.99993e-01 6.58597e-01 3.40053e-02 -3.59510e-03
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 5 ERR DEF=0.5
4.682e-03 1.464e-04 -2.241e-05 -1.666e-04 1.597e-06
1.464e-04 6.151e-03 -6.777e-04 -1.233e-04 1.335e-05
-2.241e-05 -6.777e-04 4.549e-04 1.230e-05 -4.750e-06
-1.666e-04 -1.233e-04 1.230e-05 8.594e-04 -6.227e-08
1.597e-06 1.335e-05 -4.750e-06 -6.227e-08 9.468e-06
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3 4 5
1 0.08648 1.000 0.027 -0.015 -0.083 0.008
2 0.40890 0.027 1.000 -0.405 -0.054 0.055
3 0.40825 -0.015 -0.405 1.000 0.020 -0.072
4 0.09773 -0.083 -0.054 0.020 1.000 -0.001
5 0.07801 0.008 0.055 -0.072 -0.001 1.000
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 2500
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=1949.67 FROM HESSE STATUS=OK 31 CALLS 196 TOTAL
EDM=2.47191e-05 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 a0 5.35942e-01 6.82123e-02 1.66747e-04 7.19468e-02
2 a1 1.95503e-01 7.77468e-02 2.20463e-04 -6.54793e-01
3 bkgfrac 5.40925e-01 2.12519e-02 2.38351e-04 8.19410e-02
4 mean 5.01914e+00 2.93147e-02 3.55975e-05 3.82861e-03
5 sig1frac 9.99993e-01 6.60187e-01 6.80107e-03 1.56563e+00
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 5 ERR DEF=0.5
4.682e-03 1.450e-04 -2.161e-05 -1.666e-04 -2.406e-07
1.450e-04 6.124e-03 -6.677e-04 -1.218e-04 -2.003e-06
-2.161e-05 -6.677e-04 4.519e-04 1.192e-05 7.150e-07
-1.666e-04 -1.218e-04 1.192e-05 8.594e-04 9.374e-09
-2.406e-07 -2.003e-06 7.150e-07 9.374e-09 9.576e-06
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3 4 5
1 0.08622 1.000 0.027 -0.015 -0.083 -0.001
2 0.40431 0.027 1.000 -0.401 -0.053 -0.008
3 0.40145 -0.015 -0.401 1.000 0.019 0.011
4 0.09742 -0.083 -0.053 0.019 1.000 0.000
5 0.01172 -0.001 -0.008 0.011 0.000 1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: 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
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 0.0106966
RooGaussian::sig1[ x=x mean=mean sigma=sigma1 ] = 1.19649e-05
RooGaussian::sig2[ x=x mean=mean sigma=sigma2 ] = 0.0588136
parameter snapshots
-------------------
reference_fit = (a0=0.500958 +/- 0.0231941,a1=0.160483 +/- 0.0372743,bkgfrac=0.504708 +/- 0.0113925,mean=5.01893 +/- 0.0101204,sigma1=0.5[C],sig1frac=0.818293 +/- 0.0374314,sigma2=1[C])
reference_fit_bkgonly = (a0=0.474267 +/- 0.0211215,a1=2.86676e-12 +/- 0.000176662,bkgfrac=1[C],mean=2.6195 +/- 2.39943,sigma1=0.5[C],sig1frac=0.818293 +/- 0.324473,sigma2=1[C])
named sets
----------
observables:(x)
parameters:(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2)
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.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);
void rf510_wsnamedsets()
{
// C r e a t e m o d e l a n d d a t a s e t
// -----------------------------------------------
RooWorkspace *w = new RooWorkspace("w");
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
RooDataSet *data = model->generate(*w->set("observables"), 1000);
// Fit model to data
model->fitTo(*data);
// 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
RooArgSet *params = (RooArgSet *)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
RooDataSet *refData = model.generate(x, 10000);
model.fitTo(*refData, PrintLevel(-1));
// The kTRUE 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, kTRUE);
// Make another fit with the signal component forced to zero
// and save those parameters too
bkgfrac.setVal(1);
bkgfrac.setConstant(kTRUE);
bkgfrac.removeError();
model.fitTo(*refData, PrintLevel(-1));
w.saveSnapshot("reference_fit_bkgonly", *params, kTRUE);
}
const Bool_t kTRUE
Definition RtypesCore.h:91
@ kRed
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
#define gDirectory
Definition TDirectory.h:290
#define gPad
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...
RooAbsArg * first() const
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
RooDataSet * generate(const RooArgSet &whatVars, Int_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())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:58
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.
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:121
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Chebychev polynomial p.d.f.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
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:44
TAxis * GetYaxis() const
Definition RooPlot.cxx:1263
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
The RooWorkspace is a persistable container for RooFit projects.
void Print(Option_t *opts=0) const
Print contents of the workspace.
Bool_t defineSet(const char *name, const RooArgSet &aset, Bool_t importMissing=kFALSE)
Define a named RooArgSet with given constituents.
Bool_t saveSnapshot(const char *name, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of given parameters.
Bool_t loadSnapshot(const char *name)
Load the values and attributes of the parameters in the snapshot saved with the given name.
Bool_t import(const RooAbsArg &arg, 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())
Import a RooAbsArg object, e.g.
const RooArgSet * set(const char *name)
Return pointer to previously defined named set with given nmame If no such set is found a null pointe...
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:293
The Canvas class.
Definition TCanvas.h:23
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Date
April 2009
Author
Wouter Verkerke

Definition in file rf510_wsnamedsets.C.