rf303_conditional.cxx

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #303
00004 // 
00005 // Use of tailored p.d.f as conditional p.d.fs.s
00006 // 
00007 // pdf = gauss(x,f(y),sx | y ) with f(y) = a0 + a1*y
00008 // 
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 //
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooPolyVar.h"
00021 #include "RooProdPdf.h"
00022 #include "RooPlot.h"
00023 #include "TRandom.h"
00024 #include "TCanvas.h"
00025 #include "TH1.h"
00026 using namespace RooFit ;
00027 
00028 
00029 class TestBasic303 : public RooFitTestUnit
00030 {
00031 public: 
00032 
00033 RooDataSet* makeFakeDataXY() 
00034 {
00035   RooRealVar x("x","x",-10,10) ;
00036   RooRealVar y("y","y",-10,10) ;
00037   RooArgSet coord(x,y) ;
00038 
00039   RooDataSet* d = new RooDataSet("d","d",RooArgSet(x,y)) ;
00040 
00041   for (int i=0 ; i<10000 ; i++) {
00042     Double_t tmpy = gRandom->Gaus(0,10) ;
00043     Double_t tmpx = gRandom->Gaus(0.5*tmpy,1) ;
00044     if (fabs(tmpy)<10 && fabs(tmpx)<10) {
00045       x = tmpx ;
00046       y = tmpy ;
00047       d->add(coord) ;
00048     }
00049       
00050   }
00051 
00052   return d ;
00053 }
00054 
00055 
00056 
00057   TestBasic303(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Conditional use of F(x|y)",refFile,writeRef,verbose) {} ;
00058   Bool_t testCode() {
00059 
00060   // S e t u p   c o m p o s e d   m o d e l   g a u s s ( x , m ( y ) , s )
00061   // -----------------------------------------------------------------------
00062 
00063   // Create observables
00064   RooRealVar x("x","x",-10,10) ;
00065   RooRealVar y("y","y",-10,10) ;
00066 
00067   // Create function f(y) = a0 + a1*y
00068   RooRealVar a0("a0","a0",-0.5,-5,5) ;
00069   RooRealVar a1("a1","a1",-0.5,-1,1) ;
00070   RooPolyVar fy("fy","fy",y,RooArgSet(a0,a1)) ;
00071 
00072   // Creat gauss(x,f(y),s)
00073   RooRealVar sigma("sigma","width of gaussian",0.5,0.1,2.0) ;
00074   RooGaussian model("model","Gaussian with shifting mean",x,fy,sigma) ;  
00075 
00076 
00077   // Obtain fake external experimental dataset with values for x and y
00078   RooDataSet* expDataXY = makeFakeDataXY() ;
00079 
00080 
00081 
00082   // G e n e r a t e   d a t a   f r o m   c o n d i t i o n a l   p . d . f   m o d e l ( x | y )  
00083   // ---------------------------------------------------------------------------------------------
00084 
00085   // Make subset of experimental data with only y values
00086   RooDataSet* expDataY= (RooDataSet*) expDataXY->reduce(y) ;
00087 
00088   // Generate 10000 events in x obtained from _conditional_ model(x|y) with y values taken from experimental data
00089   RooDataSet *data = model.generate(x,ProtoData(*expDataY)) ;
00090 
00091 
00092 
00093   // F i t   c o n d i t i o n a l   p . d . f   m o d e l ( x | y )   t o   d a t a
00094   // ---------------------------------------------------------------------------------------------
00095 
00096   model.fitTo(*expDataXY,ConditionalObservables(y)) ;
00097   
00098 
00099 
00100   // P r o j e c t   c o n d i t i o n a l   p . d . f   o n   x   a n d   y   d i m e n s i o n s
00101   // ---------------------------------------------------------------------------------------------
00102 
00103   // Plot x distribution of data and projection of model on x = 1/Ndata sum(data(y_i)) model(x;y_i)
00104   RooPlot* xframe = x.frame() ;
00105   expDataXY->plotOn(xframe) ;
00106   model.plotOn(xframe,ProjWData(*expDataY)) ; 
00107 
00108 
00109   // Speed up (and approximate) projection by using binned clone of data for projection
00110   RooAbsData* binnedDataY = expDataY->binnedClone() ;
00111   model.plotOn(xframe,ProjWData(*binnedDataY),LineColor(kCyan),LineStyle(kDotted),Name("Alt1")) ;
00112 
00113 
00114   // Show effect of projection with too coarse binning
00115   ((RooRealVar*)expDataY->get()->find("y"))->setBins(5) ;
00116   RooAbsData* binnedDataY2 = expDataY->binnedClone() ;
00117   model.plotOn(xframe,ProjWData(*binnedDataY2),LineColor(kRed),Name("Alt2")) ;
00118 
00119 
00120   regPlot(xframe,"rf303_plot1") ;
00121 
00122   delete binnedDataY ;
00123   delete binnedDataY2 ;
00124   delete expDataXY ;
00125   delete expDataY ;
00126   delete data ;
00127 
00128   return kTRUE ;
00129   }
00130 } ;
00131 
00132 
00133 
00134 

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