rf316_llratioplot.C

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 "RooConstVar.h"
00021 #include "RooPolynomial.h"
00022 #include "RooAddPdf.h"
00023 #include "RooProdPdf.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "RooPlot.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf316_llratioplot()
00031 {
00032 
00033   // C r e a t e   3 D   p d f   a n d   d a t a 
00034   // -------------------------------------------
00035 
00036   // Create observables
00037   RooRealVar x("x","x",-5,5) ;
00038   RooRealVar y("y","y",-5,5) ;
00039   RooRealVar z("z","z",-5,5) ;
00040 
00041   // Create signal pdf gauss(x)*gauss(y)*gauss(z) 
00042   RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
00043   RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
00044   RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
00045   RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;
00046 
00047   // Create background pdf poly(x)*poly(y)*poly(z) 
00048   RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
00049   RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
00050   RooPolynomial pz("pz","pz",z) ;
00051   RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;
00052 
00053   // Create composite pdf sig+bkg
00054   RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
00055   RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
00056 
00057   RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;
00058 
00059 
00060 
00061   // P r o j e c t   p d f   a n d   d a t a   o n   x
00062   // -------------------------------------------------
00063 
00064   // Make plain projection of data and pdf on x observable
00065   RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
00066   data->plotOn(frame) ;
00067   model.plotOn(frame) ;
00068   
00069 
00070 
00071   // 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
00072   // ----------------------------------------------------------------------------------
00073 
00074   // Calculate projection of signal and total likelihood on (y,z) observables
00075   // i.e. integrate signal and composite model over x
00076   RooAbsPdf* sigyz = sig.createProjection(x) ;
00077   RooAbsPdf* totyz = model.createProjection(x) ;
00078 
00079   // Construct the log of the signal / signal+background probability 
00080   RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;
00081 
00082 
00083 
00084   // P l o t   d a t a   w i t h   a   L L r a t i o   c u t 
00085   // -------------------------------------------------------
00086 
00087   // Calculate the llratio value for each event in the dataset
00088   data->addColumn(llratio_func) ;
00089 
00090   // Extract the subset of data with large signal likelihood
00091   RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;
00092 
00093   // Make plot frame
00094   RooPlot* frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"),Bins(40)) ;
00095 
00096   // Plot select data on frame
00097   dataSel->plotOn(frame2) ;
00098 
00099 
00100 
00101   // 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 
00102   // ---------------------------------------------------------------------------------------------
00103 
00104   // Generate large number of events for MC integration of pdf projection
00105   RooDataSet* mcprojData = model.generate(RooArgSet(x,y,z),10000) ;
00106 
00107   // Calculate LL ratio for each generated event and select MC events with llratio)0.7
00108   mcprojData->addColumn(llratio_func) ;
00109   RooDataSet* mcprojDataSel = (RooDataSet*) mcprojData->reduce(Cut("llratio>0.7")) ;
00110     
00111   // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
00112   // on set of events with the same llratio cut as was applied to data
00113   model.plotOn(frame2,ProjWData(*mcprojDataSel)) ;
00114 
00115 
00116 
00117   TCanvas* c = new TCanvas("rf316_llratioplot","rf316_llratioplot",800,400) ;
00118   c->Divide(2) ;
00119   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00120   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
00121 
00122 
00123   
00124 }

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