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

Detailed Description

View in nbviewer Open in SWAN
Organization and simultaneous fits: basic use of the 'object factory' associated with a workspace to rapidly build pdfs functions and their parameter components

#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"
using namespace RooFit;
void rf511_wsfactory_basic(bool compact = false)
{
// C r e a t i n g a n d a d d i n g b a s i c p . d . f . s
// ----------------------------------------------------------------
// Remake example pdf of tutorial rs502_wspacewrite.C:
//
// Basic pdf construction: ClassName::ObjectName(constructor arguments)
// Variable construction : VarName[x,xlo,xhi], VarName[xlo,xhi], VarName[x]
// P.d.f. addition : SUM::ObjectName(coef1*pdf1,...coefM*pdfM,pdfN)
//
if (!compact) {
// Use object factory to build pdf of tutorial rs502_wspacewrite
w->factory("Gaussian::sig1(x[-10,10],mean[5,0,10],0.5)");
w->factory("Gaussian::sig2(x,mean,1)");
w->factory("Chebychev::bkg(x,{a0[0.5,0.,1],a1[0.2,0.,1.]})");
w->factory("SUM::sig(sig1frac[0.8,0.,1.]*sig1,sig2)");
w->factory("SUM::model(bkgfrac[0.5,0.,1.]*bkg,sig)");
} else {
// Use object factory to build pdf of tutorial rs502_wspacewrite but
// - Contracted to a single line recursive expression,
// - Omitting explicit names for components that are not referred to explicitly later
w->factory("SUM::model(bkgfrac[0.5,0.,1.]*Chebychev::bkg(x[-10,10],{a0[0.5,0.,1],a1[0.2,0.,1.]}),"
"SUM(sig1frac[0.8,0.,1.]*Gaussian(x,mean[5,0,10],0.5), Gaussian(x,mean,1)))");
}
// A d v a n c e d p . d . f . c o n s t r u c t o r a r g u m e n t s
// ----------------------------------------------------------------
//
// P.d.f. constructor arguments may by any type of RooAbsArg, but also
//
// double --> converted to RooConst(...)
// {a,b,c} --> converted to RooArgSet() or RooArgList() depending on required ctor arg
// dataset name --> converted to RooAbsData reference for any dataset residing in the workspace
// enum --> Any enum label that belongs to an enum defined in the (base) class
// Make a dummy dataset pdf 'model' and import it in the workspace
RooDataSet *data = w->pdf("model")->generate(*w->var("x"), 1000);
w->import(*data, Rename("data"));
// Construct a KEYS pdf passing a dataset name and an enum type defining the
// mirroring strategy
w->factory("KeysPdf::k(x,data,NoMirror,0.2)");
// Print workspace contents
w->Print();
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
The RooWorkspace is a persistable container for RooFit projects.
RooCmdArg Rename(const char *suffix)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing dataset modelData
[#1] INFO:ObjectHandling -- RooWorkSpace::import(w) changing name of dataset from modelData to data
RooWorkspace(w) w contents
variables
---------
(a0,a1,bkgfrac,mean,sig1frac,x)
p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 0.8
RooKeysPdf::k[ x=x ] = 0.0210129
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 0.4/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 7.45331e-07/1
RooGaussian::sig1[ x=x mean=mean sigma=0.5 ] = 1.92875e-22
RooGaussian::sig2[ x=x mean=mean sigma=1 ] = 3.72665e-06
datasets
--------
RooDataSet::data(x)
Date
July 2009
Author
Wouter Verkerke

Definition in file rf511_wsfactory_basic.C.