ROOT logo
/////////////////////////////////////////////////////////////////////////
//
// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #203
// 
// Fitting and plotting in sub ranges
//
//
// 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 "RooFitResult.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit ;


void rf203_ranges()
{
  // S e t u p   m o d e l 
  // ---------------------

  // Construct observables x
  RooRealVar x("x","x",-10,10) ;
  
  // Construct gaussx(x,mx,1)
  RooRealVar mx("mx","mx",0,-10,10) ;
  RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;

  // Construct px = 1 (flat in x)
  RooPolynomial px("px","px",x) ;

  // Construct model = f*gx + (1-f)px
  RooRealVar f("f","f",0.,1.) ;
  RooAddPdf model("model","model",RooArgList(gx,px),f) ;

  // Generated 10000 events in (x,y) from p.d.f. model
  RooDataSet* modelData = model.generate(x,10000) ;

  // F i t   f u l l   r a n g e 
  // ---------------------------

  // Fit p.d.f to all data
  RooFitResult* r_full = model.fitTo(*modelData,Save(kTRUE)) ;


  // F i t   p a r t i a l   r a n g e 
  // ----------------------------------

  // Define "signal" range in x as [-3,3]
  x.setRange("signal",-3,3) ;  

  // Fit p.d.f only to data in "signal" range
  RooFitResult* r_sig = model.fitTo(*modelData,Save(kTRUE),Range("signal")) ;


  // P l o t   /   p r i n t   r e s u l t s 
  // ---------------------------------------

  // Make plot frame in x and add data and fitted model
  RooPlot* frame = x.frame(Title("Fitting a sub range")) ;
  modelData->plotOn(frame) ;
  model.plotOn(frame,Range("Full"),LineStyle(kDashed),LineColor(kRed)) ; // Add shape in full ranged dashed
  model.plotOn(frame) ; // By default only fitted range is shown
  
  // Print fit results 
  cout << "result of fit on all data " << endl ;
  r_full->Print() ;  
  cout << "result of fit in in signal region (note increased error on signal fraction)" << endl ;
  r_sig->Print() ;

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

  return ;
  
}
 rf203_ranges.C:1
 rf203_ranges.C:2
 rf203_ranges.C:3
 rf203_ranges.C:4
 rf203_ranges.C:5
 rf203_ranges.C:6
 rf203_ranges.C:7
 rf203_ranges.C:8
 rf203_ranges.C:9
 rf203_ranges.C:10
 rf203_ranges.C:11
 rf203_ranges.C:12
 rf203_ranges.C:13
 rf203_ranges.C:14
 rf203_ranges.C:15
 rf203_ranges.C:16
 rf203_ranges.C:17
 rf203_ranges.C:18
 rf203_ranges.C:19
 rf203_ranges.C:20
 rf203_ranges.C:21
 rf203_ranges.C:22
 rf203_ranges.C:23
 rf203_ranges.C:24
 rf203_ranges.C:25
 rf203_ranges.C:26
 rf203_ranges.C:27
 rf203_ranges.C:28
 rf203_ranges.C:29
 rf203_ranges.C:30
 rf203_ranges.C:31
 rf203_ranges.C:32
 rf203_ranges.C:33
 rf203_ranges.C:34
 rf203_ranges.C:35
 rf203_ranges.C:36
 rf203_ranges.C:37
 rf203_ranges.C:38
 rf203_ranges.C:39
 rf203_ranges.C:40
 rf203_ranges.C:41
 rf203_ranges.C:42
 rf203_ranges.C:43
 rf203_ranges.C:44
 rf203_ranges.C:45
 rf203_ranges.C:46
 rf203_ranges.C:47
 rf203_ranges.C:48
 rf203_ranges.C:49
 rf203_ranges.C:50
 rf203_ranges.C:51
 rf203_ranges.C:52
 rf203_ranges.C:53
 rf203_ranges.C:54
 rf203_ranges.C:55
 rf203_ranges.C:56
 rf203_ranges.C:57
 rf203_ranges.C:58
 rf203_ranges.C:59
 rf203_ranges.C:60
 rf203_ranges.C:61
 rf203_ranges.C:62
 rf203_ranges.C:63
 rf203_ranges.C:64
 rf203_ranges.C:65
 rf203_ranges.C:66
 rf203_ranges.C:67
 rf203_ranges.C:68
 rf203_ranges.C:69
 rf203_ranges.C:70
 rf203_ranges.C:71
 rf203_ranges.C:72
 rf203_ranges.C:73
 rf203_ranges.C:74
 rf203_ranges.C:75
 rf203_ranges.C:76
 rf203_ranges.C:77
 rf203_ranges.C:78
 rf203_ranges.C:79
 rf203_ranges.C:80
 rf203_ranges.C:81
 rf203_ranges.C:82
 rf203_ranges.C:83
 rf203_ranges.C:84
 rf203_ranges.C:85
 rf203_ranges.C:86
 rf203_ranges.C:87
 rf203_ranges.C:88
 rf203_ranges.C:89
 rf203_ranges.C:90