ROOT logo

From $ROOTSYS/tutorials/roofit/rf407_latextables.C

//////////////////////////////////////////////////////////////////////////
//
// 'DATA AND CATEGORIES' RooFit tutorial macro #407
// 
// Latex printing of lists and sets of RooArgSets
//
//
//
// 07/2008 - 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 "RooExponential.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf407_latextables()
{
  // S e t u p   c o m p o s i t e    p d f
  // --------------------------------------

  // Declare observable x
  RooRealVar x("x","x",0,10) ;

  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their paramaters
  RooRealVar mean("mean","mean of gaussians",5) ;
  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) ;  

  // Sum the signal components into a composite signal p.d.f.
  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
  
  // Build Chebychev polynomial p.d.f.  
  RooRealVar a0("a0","a0",0.5,0.,1.) ;
  RooRealVar a1("a1","a1",-0.2,0.,1.) ;
  RooChebychev bkg1("bkg1","Background 1",x,RooArgSet(a0,a1)) ;

  // Build expontential pdf
  RooRealVar alpha("alpha","alpha",-1) ;
  RooExponential bkg2("bkg2","Background 2",x,alpha) ;

  // Sum the background components into a composite background p.d.f.
  RooRealVar bkg1frac("sig1frac","fraction of component 1 in background",0.2,0.,1.) ;
  RooAddPdf bkg("bkg","Signal",RooArgList(bkg1,bkg2),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) ;



  // M a k e   l i s t   o f   p a r a m e t e r s   b e f o r e   a n d   a f t e r   f i t
  // ----------------------------------------------------------------------------------------

  // Make list of model parameters
  RooArgSet* params = model.getParameters(x) ;

  // Save snapshot of prefit parameters
  RooArgSet* initParams = (RooArgSet*) params->snapshot() ;

  // Do fit to data, to obtain error estimates on parameters
  RooDataSet* data = model.generate(x,1000) ;
  model.fitTo(*data) ;



  // P r i n t   l a t ex   t a b l e   o f   p a r a m e t e r s   o f   p d f 
  // --------------------------------------------------------------------------


  // Print parameter list in LaTeX for (one column with names, one column with values)
  params->printLatex() ;

  // Print parameter list in LaTeX for (names values|names values)
  params->printLatex(Columns(2)) ;

  // Print two parameter lists side by side (name values initvalues)
  params->printLatex(Sibling(*initParams)) ;

  // Print two parameter lists side by side (name values initvalues|name values initvalues)
  params->printLatex(Sibling(*initParams),Columns(2)) ;

  // Write LaTex table to file
  params->printLatex(Sibling(*initParams),OutputFile("rf407_latextables.tex")) ;


}


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