ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf705_linearmorph.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'SPECIAL PDFS' RooFit tutorial macro #705
4 //
5 // Linear interpolation between p.d.f shapes using the 'Alex Read' algorithm
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 "RooPolynomial.h"
20 #include "RooLinearMorph.h"
21 #include "RooNLLVar.h"
22 #include "TCanvas.h"
23 #include "RooPlot.h"
24 #include "TH1.h"
25 using namespace RooFit ;
26 
27 
28 // Elementary operations on a gaussian PDF
29 class TestBasic705 : public RooFitTestUnit
30 {
31 public:
32 
33  Double_t ctol() { return 5e-2 ; } // very conservative, this is a numerically difficult test
34 
35  TestBasic705(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Linear morph operator p.d.f.",refFile,writeRef,verbose) {} ;
36  Bool_t testCode() {
37 
38  // C r e a t e e n d p o i n t p d f s h a p e s
39  // ------------------------------------------------------
40 
41  // Observable
42  RooRealVar x("x","x",-20,20) ;
43 
44  // Lower end point shape: a Gaussian
45  RooRealVar g1mean("g1mean","g1mean",-10) ;
46  RooGaussian g1("g1","g1",x,g1mean,RooConst(2)) ;
47 
48  // Upper end point shape: a Polynomial
49  RooPolynomial g2("g2","g2",x,RooArgSet(RooConst(-0.03),RooConst(-0.001))) ;
50 
51 
52 
53  // C r e a t e i n t e r p o l a t i n g p d f
54  // -----------------------------------------------
55 
56  // Create interpolation variable
57  RooRealVar alpha("alpha","alpha",0,1.0) ;
58 
59  // Specify sampling density on observable and interpolation variable
60  x.setBins(1000,"cache") ;
61  alpha.setBins(50,"cache") ;
62 
63  // Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
64  // and g2(x) at a=a_max
65  RooLinearMorph lmorph("lmorph","lmorph",g1,g2,x,alpha) ;
66 
67 
68 
69  // P l o t i n t e r p o l a t i n g p d f a t v a r i o u s a l p h a
70  // -----------------------------------------------------------------------------
71 
72  // Show end points as blue curves
73  RooPlot* frame1 = x.frame() ;
74  g1.plotOn(frame1) ;
75  g2.plotOn(frame1) ;
76 
77  // Show interpolated shapes in red
78  alpha.setVal(0.125) ;
79  lmorph.plotOn(frame1,LineColor(kRed),Name("alt1")) ;
80  alpha.setVal(0.25) ;
81  lmorph.plotOn(frame1,LineColor(kRed),Name("alt2")) ;
82  alpha.setVal(0.375) ;
83  lmorph.plotOn(frame1,LineColor(kRed),Name("alt3")) ;
84  alpha.setVal(0.50) ;
85  lmorph.plotOn(frame1,LineColor(kRed),Name("alt4")) ;
86  alpha.setVal(0.625) ;
87  lmorph.plotOn(frame1,LineColor(kRed),Name("alt5")) ;
88  alpha.setVal(0.75) ;
89  lmorph.plotOn(frame1,LineColor(kRed),Name("alt6")) ;
90  alpha.setVal(0.875) ;
91  lmorph.plotOn(frame1,LineColor(kRed),Name("alt7")) ;
92  alpha.setVal(0.95) ;
93  lmorph.plotOn(frame1,LineColor(kRed),Name("alt8")) ;
94 
95 
96 
97  // S h o w 2 D d i s t r i b u t i o n o f p d f ( x , a l p h a )
98  // -----------------------------------------------------------------------
99 
100  // Create 2D histogram
101  TH1* hh = lmorph.createHistogram("hh",x,Binning(40),YVar(alpha,Binning(40))) ;
102  hh->SetLineColor(kBlue) ;
103 
104 
105  // F i t p d f t o d a t a s e t w i t h a l p h a = 0 . 8
106  // -----------------------------------------------------------------
107 
108  // Generate a toy dataset at alpha = 0.8
109  alpha=0.8 ;
110  RooDataSet* data = lmorph.generate(x,1000) ;
111 
112  // Fit pdf to toy data
113  lmorph.setCacheAlpha(kTRUE) ;
114  lmorph.fitTo(*data) ;
115 
116  // Plot fitted pdf and data overlaid
117  RooPlot* frame2 = x.frame(Bins(100)) ;
118  data->plotOn(frame2) ;
119  lmorph.plotOn(frame2) ;
120 
121 
122  // S c a n - l o g ( L ) v s a l p h a
123  // -----------------------------------------
124 
125  // Show scan -log(L) of dataset w.r.t alpha
126  RooPlot* frame3 = alpha.frame(Bins(100),Range(0.5,0.9)) ;
127 
128  // Make 2D pdf of histogram
129  RooNLLVar nll("nll","nll",lmorph,*data) ;
130  nll.plotOn(frame3,ShiftToZero()) ;
131 
132  lmorph.setCacheAlpha(kFALSE) ;
133 
134 
135  regPlot(frame1,"rf705_plot1") ;
136  regPlot(frame2,"rf705_plot2") ;
137  regPlot(frame3,"rf705_plot3") ;
138  regTH(hh,"rf705_hh") ;
139 
140  delete data ;
141 
142  return kTRUE;
143  }
144 } ;
RooCmdArg LineColor(Color_t color)
tuple g2
Definition: multifit.py:39
Definition: Rtypes.h:61
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
tuple g1
Definition: multifit.py:38
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:626
RooCmdArg ShiftToZero()
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
Class RooNLLVar implements a a -log(likelihood) calculation from a dataset and a PDF.
Definition: RooNLLVar.h:26
RooCmdArg Name(const char *name)
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
double Double_t
Definition: RtypesCore.h:55
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
The TH1 histogram class.
Definition: TH1.h:80
RooConstVar & RooConst(Double_t val)
RooCmdArg Bins(Int_t nbin)
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
Definition: Rtypes.h:61
const Bool_t kTRUE
Definition: Rtypes.h:91
RooCmdArg Binning(const RooAbsBinning &binning)