ROOT logo

From $ROOTSYS/tutorials/roofit/rf316_llratioplot.C

//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #316
// 
// Using the likelihood ratio techique to construct a signal enhanced
// one-dimensional projection of a multi-dimensional p.d.f.
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooPolynomial.h"
#include "RooAddPdf.h"
#include "RooProdPdf.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf316_llratioplot()
{

  // C r e a t e   3 D   p d f   a n d   d a t a 
  // -------------------------------------------

  // Create observables
  RooRealVar x("x","x",-5,5) ;
  RooRealVar y("y","y",-5,5) ;
  RooRealVar z("z","z",-5,5) ;

  // Create signal pdf gauss(x)*gauss(y)*gauss(z) 
  RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
  RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
  RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
  RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;

  // Create background pdf poly(x)*poly(y)*poly(z) 
  RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
  RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
  RooPolynomial pz("pz","pz",z) ;
  RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;

  // Create composite pdf sig+bkg
  RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
  RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;

  RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;



  // P r o j e c t   p d f   a n d   d a t a   o n   x
  // -------------------------------------------------

  // Make plain projection of data and pdf on x observable
  RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
  data->plotOn(frame) ;
  model.plotOn(frame) ;
  


  // 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
  // ----------------------------------------------------------------------------------

  // Calculate projection of signal and total likelihood on (y,z) observables
  // i.e. integrate signal and composite model over x
  RooAbsPdf* sigyz = sig.createProjection(x) ;
  RooAbsPdf* totyz = model.createProjection(x) ;

  // Construct the log of the signal / signal+background probability 
  RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;



  // P l o t   d a t a   w i t h   a   L L r a t i o   c u t 
  // -------------------------------------------------------

  // Calculate the llratio value for each event in the dataset
  data->addColumn(llratio_func) ;

  // Extract the subset of data with large signal likelihood
  RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;

  // Make plot frame
  RooPlot* frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"),Bins(40)) ;

  // Plot select data on frame
  dataSel->plotOn(frame2) ;



  // 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 
  // ---------------------------------------------------------------------------------------------

  // Generate large number of events for MC integration of pdf projection
  RooDataSet* mcprojData = model.generate(RooArgSet(x,y,z),10000) ;

  // Calculate LL ratio for each generated event and select MC events with llratio)0.7
  mcprojData->addColumn(llratio_func) ;
  RooDataSet* mcprojDataSel = (RooDataSet*) mcprojData->reduce(Cut("llratio>0.7")) ;
    
  // Project model on x, integrating projected observables (y,z) with Monte Carlo technique
  // on set of events with the same llratio cut as was applied to data
  model.plotOn(frame2,ProjWData(*mcprojDataSel)) ;



  TCanvas* c = new TCanvas("rf316_llratioplot","rf316_llratioplot",800,400) ;
  c->Divide(2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;


  
}
 rf316_llratioplot.C:1
 rf316_llratioplot.C:2
 rf316_llratioplot.C:3
 rf316_llratioplot.C:4
 rf316_llratioplot.C:5
 rf316_llratioplot.C:6
 rf316_llratioplot.C:7
 rf316_llratioplot.C:8
 rf316_llratioplot.C:9
 rf316_llratioplot.C:10
 rf316_llratioplot.C:11
 rf316_llratioplot.C:12
 rf316_llratioplot.C:13
 rf316_llratioplot.C:14
 rf316_llratioplot.C:15
 rf316_llratioplot.C:16
 rf316_llratioplot.C:17
 rf316_llratioplot.C:18
 rf316_llratioplot.C:19
 rf316_llratioplot.C:20
 rf316_llratioplot.C:21
 rf316_llratioplot.C:22
 rf316_llratioplot.C:23
 rf316_llratioplot.C:24
 rf316_llratioplot.C:25
 rf316_llratioplot.C:26
 rf316_llratioplot.C:27
 rf316_llratioplot.C:28
 rf316_llratioplot.C:29
 rf316_llratioplot.C:30
 rf316_llratioplot.C:31
 rf316_llratioplot.C:32
 rf316_llratioplot.C:33
 rf316_llratioplot.C:34
 rf316_llratioplot.C:35
 rf316_llratioplot.C:36
 rf316_llratioplot.C:37
 rf316_llratioplot.C:38
 rf316_llratioplot.C:39
 rf316_llratioplot.C:40
 rf316_llratioplot.C:41
 rf316_llratioplot.C:42
 rf316_llratioplot.C:43
 rf316_llratioplot.C:44
 rf316_llratioplot.C:45
 rf316_llratioplot.C:46
 rf316_llratioplot.C:47
 rf316_llratioplot.C:48
 rf316_llratioplot.C:49
 rf316_llratioplot.C:50
 rf316_llratioplot.C:51
 rf316_llratioplot.C:52
 rf316_llratioplot.C:53
 rf316_llratioplot.C:54
 rf316_llratioplot.C:55
 rf316_llratioplot.C:56
 rf316_llratioplot.C:57
 rf316_llratioplot.C:58
 rf316_llratioplot.C:59
 rf316_llratioplot.C:60
 rf316_llratioplot.C:61
 rf316_llratioplot.C:62
 rf316_llratioplot.C:63
 rf316_llratioplot.C:64
 rf316_llratioplot.C:65
 rf316_llratioplot.C:66
 rf316_llratioplot.C:67
 rf316_llratioplot.C:68
 rf316_llratioplot.C:69
 rf316_llratioplot.C:70
 rf316_llratioplot.C:71
 rf316_llratioplot.C:72
 rf316_llratioplot.C:73
 rf316_llratioplot.C:74
 rf316_llratioplot.C:75
 rf316_llratioplot.C:76
 rf316_llratioplot.C:77
 rf316_llratioplot.C:78
 rf316_llratioplot.C:79
 rf316_llratioplot.C:80
 rf316_llratioplot.C:81
 rf316_llratioplot.C:82
 rf316_llratioplot.C:83
 rf316_llratioplot.C:84
 rf316_llratioplot.C:85
 rf316_llratioplot.C:86
 rf316_llratioplot.C:87
 rf316_llratioplot.C:88
 rf316_llratioplot.C:89
 rf316_llratioplot.C:90
 rf316_llratioplot.C:91
 rf316_llratioplot.C:92
 rf316_llratioplot.C:93
 rf316_llratioplot.C:94
 rf316_llratioplot.C:95
 rf316_llratioplot.C:96
 rf316_llratioplot.C:97
 rf316_llratioplot.C:98
 rf316_llratioplot.C:99
 rf316_llratioplot.C:100
 rf316_llratioplot.C:101
 rf316_llratioplot.C:102
 rf316_llratioplot.C:103
 rf316_llratioplot.C:104
 rf316_llratioplot.C:105
 rf316_llratioplot.C:106
 rf316_llratioplot.C:107
 rf316_llratioplot.C:108
 rf316_llratioplot.C:109
 rf316_llratioplot.C:110
 rf316_llratioplot.C:111
 rf316_llratioplot.C:112
 rf316_llratioplot.C:113
 rf316_llratioplot.C:114
 rf316_llratioplot.C:115
 rf316_llratioplot.C:116
 rf316_llratioplot.C:117
 rf316_llratioplot.C:118
 rf316_llratioplot.C:119
 rf316_llratioplot.C:120
 rf316_llratioplot.C:121
 rf316_llratioplot.C:122
 rf316_llratioplot.C:123
 rf316_llratioplot.C:124
 rf316_llratioplot.C:125