rf705_linearmorph.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'SPECIAL PDFS' RooFit tutorial macro #705
00004 // 
00005 // Linear interpolation between p.d.f shapes using the 'Alex Read' algorithm
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussian.h"
00019 #include "RooPolynomial.h"
00020 #include "RooLinearMorph.h"
00021 #include "RooNLLVar.h"
00022 #include "TCanvas.h"
00023 #include "RooPlot.h"
00024 #include "TH1.h"
00025 using namespace RooFit ;
00026 
00027 
00028 // Elementary operations on a gaussian PDF
00029 class TestBasic705 : public RooFitTestUnit
00030 {
00031 public: 
00032 
00033   Double_t ctol() { return 5e-2 ; } // very conservative, this is a numerically difficult test
00034 
00035   TestBasic705(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Linear morph operator p.d.f.",refFile,writeRef,verbose) {} ;
00036   Bool_t testCode() {
00037 
00038   // C r e a t e   e n d   p o i n t   p d f   s h a p e s
00039   // ------------------------------------------------------
00040 
00041   // Observable
00042   RooRealVar x("x","x",-20,20) ;
00043 
00044   // Lower end point shape: a Gaussian
00045   RooRealVar g1mean("g1mean","g1mean",-10) ;
00046   RooGaussian g1("g1","g1",x,g1mean,RooConst(2)) ;
00047 
00048   // Upper end point shape: a Polynomial
00049   RooPolynomial g2("g2","g2",x,RooArgSet(RooConst(-0.03),RooConst(-0.001))) ;
00050 
00051 
00052  
00053   // C r e a t e   i n t e r p o l a t i n g   p d f 
00054   // -----------------------------------------------
00055 
00056   // Create interpolation variable
00057   RooRealVar alpha("alpha","alpha",0,1.0) ;
00058 
00059   // Specify sampling density on observable and interpolation variable
00060   x.setBins(1000,"cache") ;
00061   alpha.setBins(50,"cache") ;
00062 
00063   // Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
00064   // and g2(x) at a=a_max
00065   RooLinearMorph lmorph("lmorph","lmorph",g1,g2,x,alpha) ;
00066 
00067 
00068 
00069   // 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 
00070   // -----------------------------------------------------------------------------
00071 
00072   // Show end points as blue curves
00073   RooPlot* frame1 = x.frame() ;
00074   g1.plotOn(frame1) ;
00075   g2.plotOn(frame1) ;
00076 
00077   // Show interpolated shapes in red
00078   alpha.setVal(0.125) ;
00079   lmorph.plotOn(frame1,LineColor(kRed),Name("alt1")) ;
00080   alpha.setVal(0.25) ;
00081   lmorph.plotOn(frame1,LineColor(kRed),Name("alt2")) ;
00082   alpha.setVal(0.375) ;
00083   lmorph.plotOn(frame1,LineColor(kRed),Name("alt3")) ;
00084   alpha.setVal(0.50) ;
00085   lmorph.plotOn(frame1,LineColor(kRed),Name("alt4")) ;
00086   alpha.setVal(0.625) ;
00087   lmorph.plotOn(frame1,LineColor(kRed),Name("alt5")) ;
00088   alpha.setVal(0.75) ;
00089   lmorph.plotOn(frame1,LineColor(kRed),Name("alt6")) ;
00090   alpha.setVal(0.875) ;
00091   lmorph.plotOn(frame1,LineColor(kRed),Name("alt7")) ;
00092   alpha.setVal(0.95) ;
00093   lmorph.plotOn(frame1,LineColor(kRed),Name("alt8")) ;
00094 
00095 
00096 
00097   // 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 ) 
00098   // -----------------------------------------------------------------------
00099   
00100   // Create 2D histogram
00101   TH1* hh = lmorph.createHistogram("hh",x,Binning(40),YVar(alpha,Binning(40))) ;
00102   hh->SetLineColor(kBlue) ;
00103 
00104 
00105   // 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 
00106   // -----------------------------------------------------------------
00107 
00108   // Generate a toy dataset at alpha = 0.8
00109   alpha=0.8 ;
00110   RooDataSet* data = lmorph.generate(x,1000) ;
00111 
00112   // Fit pdf to toy data
00113   lmorph.setCacheAlpha(kTRUE) ;
00114   lmorph.fitTo(*data) ;
00115 
00116   // Plot fitted pdf and data overlaid
00117   RooPlot* frame2 = x.frame(Bins(100)) ;
00118   data->plotOn(frame2) ;
00119   lmorph.plotOn(frame2) ; 
00120 
00121 
00122   // S c a n   - l o g ( L )   v s   a l p h a
00123   // -----------------------------------------
00124 
00125   // Show scan -log(L) of dataset w.r.t alpha
00126   RooPlot* frame3 = alpha.frame(Bins(100),Range(0.5,0.9)) ;
00127   
00128   // Make 2D pdf of histogram  
00129   RooNLLVar nll("nll","nll",lmorph,*data) ;  
00130   nll.plotOn(frame3,ShiftToZero()) ;    
00131 
00132   lmorph.setCacheAlpha(kFALSE) ;
00133 
00134 
00135   regPlot(frame1,"rf705_plot1") ;
00136   regPlot(frame2,"rf705_plot2") ;
00137   regPlot(frame3,"rf705_plot3") ;
00138   regTH(hh,"rf705_hh") ;
00139 
00140   delete data ;
00141   
00142   return kTRUE;
00143   }
00144 } ;

Generated on Tue Jul 5 15:07:51 2011 for ROOT_528-00b_version by  doxygen 1.5.1