Logo ROOT   6.16/01
Reference Guide
rf505_asciicfg.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -nodraw
4/// 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #505
5///
6/// Reading and writing ASCII configuration files
7///
8/// \macro_output
9/// \macro_code
10/// \author 07/2008 - Wouter Verkerke
11
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooConstVar.h"
17#include "RooPolynomial.h"
18#include "RooAddPdf.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit ;
23
24
25void rf505_asciicfg()
26{
27 // C r e a t e p d f
28 // ------------------
29
30 // Construct gauss(x,m,s)
31 RooRealVar x("x","x",-10,10) ;
32 RooRealVar m("m","m",0,-10,10) ;
33 RooRealVar s("s","s",1,-10,10) ;
34 RooGaussian gauss("g","g",x,m,s) ;
35
36 // Construct poly(x,p0)
37 RooRealVar p0("p0","p0",0.01,0.,1.) ;
38 RooPolynomial poly("p","p",x,p0) ;
39
40 // Construct model = f*gauss(x) + (1-f)*poly(x)
41 RooRealVar f("f","f",0.5,0.,1.) ;
42 RooAddPdf model("model","model",RooArgSet(gauss,poly),f) ;
43
44
45
46 // F i t m o d e l t o t o y d a t a
47 // -----------------------------------------
48
49 RooDataSet* d = model.generate(x,1000) ;
50 model.fitTo(*d) ;
51
52
53 // W r i t e p a r a m e t e r s t o a s c i i f i l e
54 // -----------------------------------------------------------
55
56 // Obtain set of parameters
57 RooArgSet* params = model.getParameters(x) ;
58
59 // Write parameters to file
60 params->writeToFile("rf505_asciicfg_example.txt") ;
61
62
63 TString dir1 = gROOT->GetTutorialDir() ;
64 dir1.Append("/roofit/rf505_asciicfg.txt") ;
65 TString dir2 = gROOT->GetTutorialDir() ;
66 dir2.Append("/roofit/rf505_asciicfg_example.txt") ;
67 // R e a d p a r a m e t e r s f r o m a s c i i f i l e
68 // ----------------------------------------------------------------
69
70 // Read parameters from file
71 params->readFromFile(dir2) ;
72 params->Print("v") ;
73
74 // Read parameters from section 'Section2' of file
75 params->readFromFile(dir1,0,"Section2") ;
76 params->Print("v") ;
77
78 // Read parameters from section 'Section3' of file. Mark all
79 // variables that were processed with the "READ" attribute
80 params->readFromFile(dir1,"READ","Section3") ;
81
82 // Print the list of parameters that were not read from Section3
83 cout << "The following parameters of the were _not_ read from Section3: "
84 << (*params->selectByAttrib("READ",kFALSE)) << endl ;
85
86
87 // Read parameters from section 'Section4' of file, which contains
88 // 'include file' statement of rf505_asciicfg_example.txt
89 // so that we effective read the same
90
91 params->readFromFile(dir1,0,"Section4") ;
92 params->Print("v") ;
93
94
95
96}
#define d(i)
Definition: RSha256.hxx:102
#define f(i)
Definition: RSha256.hxx:104
const Bool_t kFALSE
Definition: RtypesCore.h:88
#define gROOT
Definition: TROOT.h:410
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
RooAbsCollection * selectByAttrib(const char *name, Bool_t value) const
Create a subset of the current collection, consisting only of those elements with the specified attri...
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Bool_t readFromFile(const char *fileName, const char *flagReadAtt=0, const char *section=0, Bool_t verbose=kFALSE)
Read contents of the argset from specified file.
Definition: RooArgSet.cxx:660
void writeToFile(const char *fileName) const
Write contents of the argset to specified file.
Definition: RooArgSet.cxx:644
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooPolynomial implements a polynomial p.d.f of the form.
Definition: RooPolynomial.h:28
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
Double_t x[n]
Definition: legend1.C:17
static constexpr double s
static constexpr double gauss
auto * m
Definition: textangle.C:8