rf903_numintcache.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'NUMERIC ALGORITHM TUNING' RooFit tutorial macro #903 
00004 // 
00005 //  Caching of slow numeric integrals and parameterizations of slow
00006 //  numeric integrals
00007 //
00008 // 07/2008 - Wouter Verkerke 
00009 // 
00010 /////////////////////////////////////////////////////////////////////////
00011 
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooDataHist.h"
00018 #include "RooGaussian.h"
00019 #include "RooConstVar.h"
00020 #include "TCanvas.h"
00021 #include "TAxis.h"
00022 #include "RooPlot.h"
00023 #include "RooWorkspace.h"
00024 #include "RooExpensiveObjectCache.h"
00025 #include "TFile.h"
00026 #include "TH1.h"
00027 
00028 using namespace RooFit ;
00029 
00030 RooWorkspace* getWorkspace(Int_t mode) ;
00031 
00032 void rf903_numintcache(Int_t mode=0)
00033 {
00034   // Mode = 0 : Run plain fit (slow)
00035   // Mode = 1 : Generate workspace with precalculated integral and store it on file (prepare for accelerated running)
00036   // Mode = 2 : Run fit from previously stored workspace including cached integrals (fast, requires run in mode=1 first)
00037 
00038   // C r e a t e ,   s a v e   o r   l o a d   w o r k s p a c e   w i t h   p . d . f . 
00039   // -----------------------------------------------------------------------------------
00040 
00041   // Make/load workspace, exit here in mode 1
00042   RooWorkspace* w = getWorkspace(mode) ;
00043   if (mode==1) {
00044 
00045     // Show workspace that was created
00046     w->Print() ;
00047 
00048     // Show plot of cached integral values
00049     RooDataHist* hhcache = (RooDataHist*) w->expensiveObjectCache().getObj(1) ;
00050 
00051     new TCanvas("rf903_numintcache","rf903_numintcache",600,600) ;
00052     hhcache->createHistogram("a")->Draw() ;
00053     
00054     return ;
00055   }
00056 
00057 
00058   // U s e   p . d . f .   f r o m   w o r k s p a c e   f o r   g e n e r a t i o n   a n d   f i t t i n g 
00059   // -----------------------------------------------------------------------------------
00060 
00061   // This is always slow (need to find maximum function value empirically in 3D space)
00062   RooDataSet* d = w->pdf("model")->generate(RooArgSet(*w->var("x"),*w->var("y"),*w->var("z")),1000) ;
00063 
00064   // This is slow in mode 0, but fast in mode 1
00065   w->pdf("model")->fitTo(*d,Verbose(kTRUE),Timer(kTRUE)) ; 
00066 
00067   // Projection on x (always slow as 2D integral over Y,Z at fitted value of a is not cached)
00068   RooPlot* framex = w->var("x")->frame(Title("Projection of 3D model on X")) ;
00069   d->plotOn(framex) ;
00070   w->pdf("model")->plotOn(framex) ;
00071 
00072   // Draw x projection on canvas
00073   new TCanvas("rf903_numintcache","rf903_numintcache",600,600) ;
00074   framex->Draw() ;
00075 
00076   // Make workspace available on command line after macro finishes
00077   gDirectory->Add(w) ;
00078 
00079   return ;
00080 
00081  
00082 }
00083 
00084 
00085 
00086 RooWorkspace* getWorkspace(Int_t mode) 
00087 {
00088   // C r e a t e ,   s a v e   o r   l o a d   w o r k s p a c e   w i t h   p . d . f . 
00089   // -----------------------------------------------------------------------------------
00090   //
00091   // Mode = 0 : Create workspace for plain running (no integral caching)
00092   // Mode = 1 : Generate workspace with precalculated integral and store it on file
00093   // Mode = 2 : Load previously stored workspace from file
00094 
00095   RooWorkspace* w(0) ;
00096 
00097   if (mode!=2) {
00098 
00099     // Create empty workspace workspace 
00100     w = new RooWorkspace("w",1) ;
00101 
00102     // Make a difficult to normalize  p.d.f. in 3 dimensions that is integrated numerically.
00103     w->factory("EXPR::model('1/((x-a)*(x-a)+0.01)+1/((y-a)*(y-a)+0.01)+1/((z-a)*(z-a)+0.01)',x[-1,1],y[-1,1],z[-1,1],a[-5,5])") ;
00104   }
00105 
00106   if (mode==1) {
00107     
00108     // Instruct model to precalculate normalization integral that integrate at least
00109     // two dimensions numerically. In this specific case the integral value for
00110     // all values of parameter 'a' are stored in a histogram and available for use 
00111     // in subsequent fitting and plotting operations (interpolation is applied)
00112     w->pdf("model")->setNormValueCaching(3) ;
00113     
00114     // Evaluate p.d.f. once to trigger filling of cache
00115     RooArgSet normSet(*w->var("x"),*w->var("y"),*w->var("z")) ;
00116     w->pdf("model")->getVal(&normSet) ;
00117     w->writeToFile("rf903_numintcache.root") ;
00118 
00119   } 
00120 
00121   if (mode==2) {    
00122     // Load preexisting workspace from file in mode==2
00123     TFile* f = new TFile("rf903_numintcache.root") ;
00124     w = (RooWorkspace*) f->Get("w") ;
00125   }
00126 
00127   // Return created or loaded workspace
00128   return w ;
00129 }

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