Logo ROOT   6.10/09
Reference Guide
rf607_fitresult.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #607
4 //
5 // Demonstration of options of the RooFitResult class
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooAddPdf.h"
20 #include "RooChebychev.h"
21 #include "RooFitResult.h"
22 #include "TCanvas.h"
23 #include "RooPlot.h"
24 #include "TFile.h"
25 #include "TStyle.h"
26 #include "TH2.h"
27 
28 using namespace RooFit ;
29 
30 
31 class TestBasic607 : public RooFitTestUnit
32 {
33 public:
34  TestBasic607(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fit Result functionality",refFile,writeRef,verbose) {} ;
35  Bool_t testCode() {
36 
37  // C r e a t e p d f , d a t a
38  // --------------------------------
39 
40  // Declare observable x
41  RooRealVar x("x","x",0,10) ;
42 
43  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
44  RooRealVar mean("mean","mean of gaussians",5,-10,10) ;
45  RooRealVar sigma1("sigma1","width of gaussians",0.5,0.1,10) ;
46  RooRealVar sigma2("sigma2","width of gaussians",1,0.1,10) ;
47 
48  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
49  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
50 
51  // Build Chebychev polynomial p.d.f.
52  RooRealVar a0("a0","a0",0.5,0.,1.) ;
53  RooRealVar a1("a1","a1",-0.2) ;
54  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
55 
56  // Sum the signal components into a composite signal p.d.f.
57  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
58  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
59 
60  // Sum the composite signal and background
61  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
62  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
63 
64  // Generate 1000 events
65  RooDataSet* data = model.generate(x,1000) ;
66 
67 
68 
69  // F i t p d f t o d a t a , s a v e f i t r e s u l t
70  // -------------------------------------------------------------
71 
72  // Perform fit and save result
73  RooFitResult* r = model.fitTo(*data,Save()) ;
74 
75 
76  // V i s u a l i z e c o r r e l a t i o n m a t r i x
77  // -------------------------------------------------------
78 
79  // Construct 2D color plot of correlation matrix
80  gStyle->SetOptStat(0) ;
81  gStyle->SetPalette(1) ;
82  TH2* hcorr = r->correlationHist() ;
83 
84 
85  // Sample dataset with parameter values according to distribution
86  // of covariance matrix of fit result
87  RooDataSet randPars("randPars","randPars",r->floatParsFinal()) ;
88  for (Int_t i=0 ; i<10000 ; i++) {
89  randPars.add(r->randomizePars()) ;
90  }
91 
92  // make histogram of 2D distribution in sigma1 vs sig1frac
93  TH1* hhrand = randPars.createHistogram("hhrand",sigma1,Binning(35,0.25,0.65),YVar(sig1frac,Binning(35,0.3,1.1))) ;
94 
95  regTH(hcorr,"rf607_hcorr") ;
96  regTH(hhrand,"rf607_hhand") ;
97 
98  delete data ;
99  delete r ;
100 
101  return kTRUE ;
102 
103  }
104 } ;
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
const RooArgList & floatParsFinal() const
Definition: RooFitResult.h:108
const RooArgList & randomizePars() const
Return a list of floating parameter values that are perturbed from the final fit values by random amo...
R__EXTERN TStyle * gStyle
Definition: TStyle.h:402
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
TH2 * correlationHist(const char *name="correlation_matrix") const
Return TH2D of correlation matrix.
TRandom2 r(17)
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0)
Add a data point, with its coordinates specified in the &#39;data&#39; argset, to the data set...
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
The TH1 histogram class.
Definition: TH1.h:56
RooCmdArg Save(Bool_t flag=kTRUE)
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1267
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const Bool_t kTRUE
Definition: RtypesCore.h:91
void SetPalette(Int_t ncolors=kBird, Int_t *colors=0, Float_t alpha=1.)
See TColor::SetPalette.
Definition: TStyle.cxx:1460
RooCmdArg Binning(const RooAbsBinning &binning)