rf705_linearmorph.C

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 "RooConstVar.h"
00020 #include "RooPolynomial.h"
00021 #include "RooIntegralMorph.h"
00022 #include "RooNLLVar.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 #include "TH1.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf705_linearmorph()
00031 {
00032   // C r e a t e   e n d   p o i n t   p d f   s h a p e s
00033   // ------------------------------------------------------
00034 
00035   // Observable
00036   RooRealVar x("x","x",-20,20) ;
00037 
00038   // Lower end point shape: a Gaussian
00039   RooRealVar g1mean("g1mean","g1mean",-10) ;
00040   RooGaussian g1("g1","g1",x,g1mean,RooConst(2)) ;
00041 
00042   // Upper end point shape: a Polynomial
00043   RooPolynomial g2("g2","g2",x,RooArgSet(RooConst(-0.03),RooConst(-0.001))) ;
00044 
00045 
00046  
00047   // C r e a t e   i n t e r p o l a t i n g   p d f 
00048   // -----------------------------------------------
00049 
00050   // Create interpolation variable
00051   RooRealVar alpha("alpha","alpha",0,1.0) ;
00052 
00053   // Specify sampling density on observable and interpolation variable
00054   x.setBins(1000,"cache") ;
00055   alpha.setBins(50,"cache") ;
00056 
00057   // Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
00058   // and g2(x) at a=a_max
00059   RooIntegralMorph lmorph("lmorph","lmorph",g1,g2,x,alpha) ;
00060 
00061 
00062 
00063   // 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 
00064   // -----------------------------------------------------------------------------
00065 
00066   // Show end points as blue curves
00067   RooPlot* frame1 = x.frame() ;
00068   g1.plotOn(frame1) ;
00069   g2.plotOn(frame1) ;
00070 
00071   // Show interpolated shapes in red
00072   alpha.setVal(0.125) ;
00073   lmorph.plotOn(frame1,LineColor(kRed)) ;
00074   alpha.setVal(0.25) ;
00075   lmorph.plotOn(frame1,LineColor(kRed)) ;
00076   alpha.setVal(0.375) ;
00077   lmorph.plotOn(frame1,LineColor(kRed)) ;
00078   alpha.setVal(0.50) ;
00079   lmorph.plotOn(frame1,LineColor(kRed)) ;
00080   alpha.setVal(0.625) ;
00081   lmorph.plotOn(frame1,LineColor(kRed)) ;
00082   alpha.setVal(0.75) ;
00083   lmorph.plotOn(frame1,LineColor(kRed)) ;
00084   alpha.setVal(0.875) ;
00085   lmorph.plotOn(frame1,LineColor(kRed)) ;
00086   alpha.setVal(0.95) ;
00087   lmorph.plotOn(frame1,LineColor(kRed)) ;
00088 
00089 
00090 
00091   // 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 ) 
00092   // -----------------------------------------------------------------------
00093   
00094   // Create 2D histogram
00095   TH1* hh = lmorph.createHistogram("hh",x,Binning(40),YVar(alpha,Binning(40))) ;
00096   hh->SetLineColor(kBlue) ;
00097 
00098 
00099   // 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 
00100   // -----------------------------------------------------------------
00101 
00102   // Generate a toy dataset at alpha = 0.8
00103   alpha=0.8 ;
00104   RooDataSet* data = lmorph.generate(x,1000) ;
00105 
00106   // Fit pdf to toy data
00107   lmorph.setCacheAlpha(kTRUE) ;
00108   lmorph.fitTo(*data,Verbose(kTRUE)) ;
00109 
00110   // Plot fitted pdf and data overlaid
00111   RooPlot* frame2 = x.frame(Bins(100)) ;
00112   data->plotOn(frame2) ;
00113   lmorph.plotOn(frame2) ; 
00114 
00115 
00116   // S c a n   - l o g ( L )   v s   a l p h a
00117   // -----------------------------------------
00118 
00119   // Show scan -log(L) of dataset w.r.t alpha
00120   RooPlot* frame3 = alpha.frame(Bins(100),Range(0.1,0.9)) ;
00121   
00122   // Make 2D pdf of histogram  
00123   RooNLLVar nll("nll","nll",lmorph,*data) ;  
00124   nll.plotOn(frame3,ShiftToZero()) ;    
00125 
00126   lmorph.setCacheAlpha(kFALSE) ;
00127 
00128 
00129 
00130   TCanvas* c = new TCanvas("rf705_linearmorph","rf705_linearmorph",800,800) ;
00131   c->Divide(2,2) ;
00132   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.6) ; frame1->Draw() ;
00133   c->cd(2) ; gPad->SetLeftMargin(0.20) ; hh->GetZaxis()->SetTitleOffset(2.5) ; hh->Draw("surf") ;
00134   c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.4) ; frame3->Draw() ;
00135   c->cd(4) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
00136   
00137   
00138   return ;
00139 
00140 }

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