ROOT logo

From $ROOTSYS/tutorials/roofit/rf511_wsfactory_basic.C

/////////////////////////////////////////////////////////////////////////
//
// 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #511
// 
//  Basic use of the 'object factory' associated with a workspace
//  to rapidly build p.d.f.s functions and their parameter components
//
//
// 04/2009 - Wouter Verkerke 
//
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#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"
using namespace RooFit ;


void rf511_wsfactory_basic(Bool_t compact=kFALSE)
{
  RooWorkspace* w = new RooWorkspace("w") ;

  // 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 p.d.f. of tutorial rs502_wspacewrite.C:
  //
  // Basic p.d.f. 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 p.d.f. 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 p.d.f. 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_t --> converted to RooConst(...)
  // {a,b,c} --> converted to RooArgSet() or RooArgList() depending on required ctor arg
  // dataset name --> convered 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 p.d.f. '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 p.d.f. 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() ;

  // Make workspace visible on command line
  gDirectory->Add(w) ;

}


 rf511_wsfactory_basic.C:1
 rf511_wsfactory_basic.C:2
 rf511_wsfactory_basic.C:3
 rf511_wsfactory_basic.C:4
 rf511_wsfactory_basic.C:5
 rf511_wsfactory_basic.C:6
 rf511_wsfactory_basic.C:7
 rf511_wsfactory_basic.C:8
 rf511_wsfactory_basic.C:9
 rf511_wsfactory_basic.C:10
 rf511_wsfactory_basic.C:11
 rf511_wsfactory_basic.C:12
 rf511_wsfactory_basic.C:13
 rf511_wsfactory_basic.C:14
 rf511_wsfactory_basic.C:15
 rf511_wsfactory_basic.C:16
 rf511_wsfactory_basic.C:17
 rf511_wsfactory_basic.C:18
 rf511_wsfactory_basic.C:19
 rf511_wsfactory_basic.C:20
 rf511_wsfactory_basic.C:21
 rf511_wsfactory_basic.C:22
 rf511_wsfactory_basic.C:23
 rf511_wsfactory_basic.C:24
 rf511_wsfactory_basic.C:25
 rf511_wsfactory_basic.C:26
 rf511_wsfactory_basic.C:27
 rf511_wsfactory_basic.C:28
 rf511_wsfactory_basic.C:29
 rf511_wsfactory_basic.C:30
 rf511_wsfactory_basic.C:31
 rf511_wsfactory_basic.C:32
 rf511_wsfactory_basic.C:33
 rf511_wsfactory_basic.C:34
 rf511_wsfactory_basic.C:35
 rf511_wsfactory_basic.C:36
 rf511_wsfactory_basic.C:37
 rf511_wsfactory_basic.C:38
 rf511_wsfactory_basic.C:39
 rf511_wsfactory_basic.C:40
 rf511_wsfactory_basic.C:41
 rf511_wsfactory_basic.C:42
 rf511_wsfactory_basic.C:43
 rf511_wsfactory_basic.C:44
 rf511_wsfactory_basic.C:45
 rf511_wsfactory_basic.C:46
 rf511_wsfactory_basic.C:47
 rf511_wsfactory_basic.C:48
 rf511_wsfactory_basic.C:49
 rf511_wsfactory_basic.C:50
 rf511_wsfactory_basic.C:51
 rf511_wsfactory_basic.C:52
 rf511_wsfactory_basic.C:53
 rf511_wsfactory_basic.C:54
 rf511_wsfactory_basic.C:55
 rf511_wsfactory_basic.C:56
 rf511_wsfactory_basic.C:57
 rf511_wsfactory_basic.C:58
 rf511_wsfactory_basic.C:59
 rf511_wsfactory_basic.C:60
 rf511_wsfactory_basic.C:61
 rf511_wsfactory_basic.C:62
 rf511_wsfactory_basic.C:63
 rf511_wsfactory_basic.C:64
 rf511_wsfactory_basic.C:65
 rf511_wsfactory_basic.C:66
 rf511_wsfactory_basic.C:67
 rf511_wsfactory_basic.C:68
 rf511_wsfactory_basic.C:69
 rf511_wsfactory_basic.C:70
 rf511_wsfactory_basic.C:71
 rf511_wsfactory_basic.C:72
 rf511_wsfactory_basic.C:73
 rf511_wsfactory_basic.C:74
 rf511_wsfactory_basic.C:75
 rf511_wsfactory_basic.C:76
 rf511_wsfactory_basic.C:77
 rf511_wsfactory_basic.C:78
 rf511_wsfactory_basic.C:79
 rf511_wsfactory_basic.C:80
 rf511_wsfactory_basic.C:81
 rf511_wsfactory_basic.C:82
 rf511_wsfactory_basic.C:83
 rf511_wsfactory_basic.C:84
 rf511_wsfactory_basic.C:85
 rf511_wsfactory_basic.C:86
 rf511_wsfactory_basic.C:87
 rf511_wsfactory_basic.C:88
 rf511_wsfactory_basic.C:89
 rf511_wsfactory_basic.C:90
 rf511_wsfactory_basic.C:91
 rf511_wsfactory_basic.C:92
 rf511_wsfactory_basic.C:93