Logo ROOT   6.07/09
Reference Guide
rf407_latextables.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'DATA AND CATEGORIES' RooFit tutorial macro #407
5 ///
6 /// Latex printing of lists and sets of RooArgSets
7 ///
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 #include "RooRealVar.h"
13 #include "RooDataSet.h"
14 #include "RooGaussian.h"
15 #include "RooConstVar.h"
16 #include "RooChebychev.h"
17 #include "RooAddPdf.h"
18 #include "RooExponential.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "RooPlot.h"
22 using namespace RooFit ;
23 
24 
25 void rf407_latextables()
26 {
27  // S e t u p c o m p o s i t e p d f
28  // --------------------------------------
29 
30  // Declare observable x
31  RooRealVar x("x","x",0,10) ;
32 
33  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
34  RooRealVar mean("mean","mean of gaussians",5) ;
35  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
36  RooRealVar sigma2("sigma2","width of gaussians",1) ;
37  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
38  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
39 
40  // Sum the signal components into a composite signal p.d.f.
41  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
42  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
43 
44  // Build Chebychev polynomial p.d.f.
45  RooRealVar a0("a0","a0",0.5,0.,1.) ;
46  RooRealVar a1("a1","a1",0.2,0.,1.) ;
47  RooChebychev bkg1("bkg1","Background 1",x,RooArgSet(a0,a1)) ;
48 
49  // Build expontential pdf
50  RooRealVar alpha("alpha","alpha",-1) ;
51  RooExponential bkg2("bkg2","Background 2",x,alpha) ;
52 
53  // Sum the background components into a composite background p.d.f.
54  RooRealVar bkg1frac("sig1frac","fraction of component 1 in background",0.2,0.,1.) ;
55  RooAddPdf bkg("bkg","Signal",RooArgList(bkg1,bkg2),sig1frac) ;
56 
57  // Sum the composite signal and background
58  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
59  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
60 
61 
62 
63  // 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
64  // ----------------------------------------------------------------------------------------
65 
66  // Make list of model parameters
67  RooArgSet* params = model.getParameters(x) ;
68 
69  // Save snapshot of prefit parameters
70  RooArgSet* initParams = (RooArgSet*) params->snapshot() ;
71 
72  // Do fit to data, to obtain error estimates on parameters
73  RooDataSet* data = model.generate(x,1000) ;
74  model.fitTo(*data) ;
75 
76 
77 
78  // 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
79  // --------------------------------------------------------------------------
80 
81 
82  // Print parameter list in LaTeX for (one column with names, one column with values)
83  params->printLatex() ;
84 
85  // Print parameter list in LaTeX for (names values|names values)
86  params->printLatex(Columns(2)) ;
87 
88  // Print two parameter lists side by side (name values initvalues)
89  params->printLatex(Sibling(*initParams)) ;
90 
91  // Print two parameter lists side by side (name values initvalues|name values initvalues)
92  params->printLatex(Sibling(*initParams),Columns(2)) ;
93 
94  // Write LaTex table to file
95  params->printLatex(Sibling(*initParams),OutputFile("rf407_latextables.tex")) ;
96 
97 
98 }
99 
100 
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
Exponential p.d.f.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooCmdArg Columns(Int_t ncol)
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooCmdArg Sibling(const RooAbsCollection &sibling)
RooCmdArg OutputFile(const char *fileName)
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
void printLatex(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
Output content of collection as LaTex table.