rf316_llratioplot.cxx

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #316
00004 // 
00005 // Using the likelihood ratio techique to construct a signal enhanced
00006 // one-dimensional projection of a multi-dimensional p.d.f.
00007 //
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 "RooPolynomial.h"
00021 #include "RooAddPdf.h"
00022 #include "RooProdPdf.h"
00023 #include "TCanvas.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026 
00027 
00028 class TestBasic316 : public RooFitTestUnit
00029 {
00030 public: 
00031   TestBasic316(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Likelihood ratio projection plot",refFile,writeRef,verbose) {} ;
00032   Bool_t testCode() {
00033 
00034   // C r e a t e   3 D   p d f   a n d   d a t a 
00035   // -------------------------------------------
00036 
00037   // Create observables
00038   RooRealVar x("x","x",-5,5) ;
00039   RooRealVar y("y","y",-5,5) ;
00040   RooRealVar z("z","z",-5,5) ;
00041 
00042   // Create signal pdf gauss(x)*gauss(y)*gauss(z) 
00043   RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
00044   RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
00045   RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
00046   RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;
00047 
00048   // Create background pdf poly(x)*poly(y)*poly(z) 
00049   RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
00050   RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
00051   RooPolynomial pz("pz","pz",z) ;
00052   RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;
00053 
00054   // Create composite pdf sig+bkg
00055   RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
00056   RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
00057 
00058   RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;
00059 
00060 
00061 
00062   // P r o j e c t   p d f   a n d   d a t a   o n   x
00063   // -------------------------------------------------
00064 
00065   // Make plain projection of data and pdf on x observable
00066   RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
00067   data->plotOn(frame) ;
00068   model.plotOn(frame) ;
00069   
00070 
00071 
00072   // D e f i n e   p r o j e c t e d   s i g n a l   l i k e l i h o o d   r a t i o
00073   // ----------------------------------------------------------------------------------
00074 
00075   // Calculate projection of signal and total likelihood on (y,z) observables
00076   // i.e. integrate signal and composite model over x
00077   RooAbsPdf* sigyz = sig.createProjection(x) ;
00078   RooAbsPdf* totyz = model.createProjection(x) ;
00079 
00080   // Construct the log of the signal / signal+background probability 
00081   RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;
00082 
00083 
00084 
00085   // P l o t   d a t a   w i t h   a   L L r a t i o   c u t 
00086   // -------------------------------------------------------
00087 
00088   // Calculate the llratio value for each event in the dataset
00089   data->addColumn(llratio_func) ;
00090 
00091   // Extract the subset of data with large signal likelihood
00092   RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;
00093 
00094   // Make plot frame
00095   RooPlot* frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"),Bins(40)) ;
00096 
00097   // Plot select data on frame
00098   dataSel->plotOn(frame2) ;
00099 
00100 
00101 
00102   // M a k e   M C   p r o j e c t i o n   o f   p d f   w i t h   s a m e   L L r a t i o   c u t 
00103   // ---------------------------------------------------------------------------------------------
00104 
00105   // Generate large number of events for MC integration of pdf projection
00106   RooDataSet* mcprojData = model.generate(RooArgSet(x,y,z),10000) ;
00107 
00108   // Calculate LL ratio for each generated event and select MC events with llratio)0.7
00109   mcprojData->addColumn(llratio_func) ;
00110   RooDataSet* mcprojDataSel = (RooDataSet*) mcprojData->reduce(Cut("llratio>0.7")) ;
00111     
00112   // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
00113   // on set of events with the same llratio cut as was applied to data
00114   model.plotOn(frame2,ProjWData(*mcprojDataSel)) ;
00115 
00116 
00117   regPlot(frame,"rf316_plot1") ;
00118   regPlot(frame2,"rf316_plot2") ;
00119   
00120   delete data ;
00121   delete dataSel ;
00122   delete mcprojData ;
00123   delete mcprojDataSel ;
00124   delete sigyz ;
00125   delete totyz ;
00126 
00127   return kTRUE ;
00128   }  
00129 } ;

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