From $ROOTSYS/tutorials/roofit/rf202_extendedmlfit.C

//////////////////////////////////////////////////////////////////////////
//
// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #202
// 
// Setting up an extended maximum likelihood fit
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

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


void rf202_extendedmlfit()
{

  // S e t u p   c o m p o n e n t   p d f s 
  // ---------------------------------------

  // Declare observable x
  RooRealVar x("x","x",0,10) ;

  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
  RooRealVar mean("mean","mean of gaussians",5) ;
  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
  RooRealVar sigma2("sigma2","width of gaussians",1) ;

  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;  
  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;  
  
  // Build Chebychev polynomial p.d.f.  
  RooRealVar a0("a0","a0",0.5,0.,1.) ;
  RooRealVar a1("a1","a1",-0.2,0.,1.) ;
  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;

  // Sum the signal components into a composite signal p.d.f.
  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;

  /////////////////////
  // M E T H O D   1 //
  /////////////////////


  // C o n s t r u c t   e x t e n d e d   c o m p o s i t e   m o d e l 
  // -------------------------------------------------------------------

  // Sum the composite signal and background into an extended pdf nsig*sig+nbkg*bkg
  RooRealVar nsig("nsig","number of signal events",500,0.,10000) ;
  RooRealVar nbkg("nbkg","number of background events",500,0,10000) ;
  RooAddPdf  model("model","(g1+g2)+a",RooArgList(bkg,sig),RooArgList(nbkg,nsig)) ;



  // S a m p l e ,   f i t   a n d   p l o t   e x t e n d e d   m o d e l 
  // ---------------------------------------------------------------------

  // Generate a data sample of expected number events in x from model
  // = model.expectedEvents() = nsig+nbkg
  RooDataSet *data = model.generate(x) ;

  // Fit model to data, extended ML term automatically included
  model.fitTo(*data) ; 
  
  // Plot data and PDF overlaid, use expected number of events for p.d.f projection normalization
  // rather than observed number of events (==data->numEntries())
  RooPlot* xframe = x.frame(Title("extended ML fit example")) ;
  data->plotOn(xframe) ;
  model.plotOn(xframe,Normalization(1.0,RooAbsReal::RelativeExpected)) ;

  // Overlay the background component of model with a dashed line
  model.plotOn(xframe,Components(bkg),LineStyle(kDashed),Normalization(1.0,RooAbsReal::RelativeExpected)) ;

  // Overlay the background+sig2 components of model with a dotted line
  model.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineStyle(kDotted),Normalization(1.0,RooAbsReal::RelativeExpected)) ;

  // Print structure of composite p.d.f.
  model.Print("t") ;


  /////////////////////
  // M E T H O D   2 //
  /////////////////////

  // C o n s t r u c t   e x t e n d e d   c o m p o n e n t s   f i r s t
  // ---------------------------------------------------------------------

  // Associated nsig/nbkg as expected number of events with sig/bkg
  RooExtendPdf esig("esig","extended signal p.d.f",sig,nsig) ;
  RooExtendPdf ebkg("ebkg","extended background p.d.f",bkg,nbkg) ;


  // S u m   e x t e n d e d   c o m p o n e n t s   w i t h o u t   c o e f s 
  // -------------------------------------------------------------------------

  // Construct sum of two extended p.d.f. (no coefficients required)
  RooAddPdf  model2("model2","(g1+g2)+a",RooArgList(ebkg,esig)) ;



  // Draw the frame on the canvas
  new TCanvas("rf202_composite","rf202_composite",600,600) ;
  gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;

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