ROOT logo

From $ROOTSYS/tutorials/roofit/rf314_paramfitrange.C

//////////////////////////////////////////////////////////////////////////
//
// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #314
// 
// Working with parameterized ranges in a fit. This an example of a
// fit with an acceptance that changes per-event 
//
//  pdf = exp(-t/tau) with t[tmin,5]
//
//  where t and tmin are both observables in the dataset
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooExponential.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "RooFitResult.h"

using namespace RooFit ;


void rf314_paramfitrange()
{

  // D e f i n e   o b s e r v a b l e s   a n d   d e c a y   p d f 
  // ---------------------------------------------------------------

  // Declare observables
  RooRealVar t("t","t",0,5) ;
  RooRealVar tmin("tmin","tmin",0,0,5) ;

  // Make parameterized range in t : [tmin,5]
  t.setRange(tmin,RooConst(t.getMax())) ;

  // Make pdf
  RooRealVar tau("tau","tau",-1.54,-10,-0.1) ;
  RooExponential model("model","model",t,tau) ;



  // C r e a t e   i n p u t   d a t a 
  // ------------------------------------

  // Generate complete dataset without acceptance cuts (for reference)
  RooDataSet* dall = model.generate(t,10000) ;

  // Generate a (fake) prototype dataset for acceptance limit values
  RooDataSet* tmp = RooGaussian("gmin","gmin",tmin,RooConst(0),RooConst(0.5)).generate(tmin,5000) ;

  // Generate dataset with t values that observe (t>tmin)
  RooDataSet* dacc = model.generate(t,ProtoData(*tmp)) ;



  // F i t   p d f   t o   d a t a   i n   a c c e p t a n c e   r e g i o n
  // -----------------------------------------------------------------------

  RooFitResult* r = model.fitTo(*dacc,Save()) ;


 
  // P l o t   f i t t e d   p d f   o n   f u l l   a n d   a c c e p t e d   d a t a 
  // ---------------------------------------------------------------------------------

  // Make plot frame, add datasets and overlay model
  RooPlot* frame = t.frame(Title("Fit to data with per-event acceptance")) ;
  dall->plotOn(frame,MarkerColor(kRed),LineColor(kRed)) ;
  model.plotOn(frame) ;
  dacc->plotOn(frame) ;

  // Print fit results to demonstrate absence of bias
  r->Print("v") ;


  new TCanvas("rf314_paramranges","rf314_paramranges",600,600) ;
  gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;

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