Logo ROOT   6.16/01
Reference Guide
rf314_paramfitrange.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #314
5///
6/// Working with parametrized ranges in a fit. This an example of a
7/// fit with an acceptance that changes per-event
8///
9/// pdf = exp(-t/tau) with t[tmin,5]
10///
11/// where t and tmin are both observables in the dataset
12///
13/// \macro_image
14/// \macro_output
15/// \macro_code
16/// \author 07/2008 - Wouter Verkerke
17
18
19#include "RooRealVar.h"
20#include "RooDataSet.h"
21#include "RooGaussian.h"
22#include "RooConstVar.h"
23#include "RooExponential.h"
24#include "TCanvas.h"
25#include "TAxis.h"
26#include "RooPlot.h"
27#include "RooFitResult.h"
28
29using namespace RooFit ;
30
31
33{
34
35 // 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
36 // ---------------------------------------------------------------
37
38 // Declare observables
39 RooRealVar t("t","t",0,5) ;
40 RooRealVar tmin("tmin","tmin",0,0,5) ;
41
42 // Make parametrized range in t : [tmin,5]
43 t.setRange(tmin,RooConst(t.getMax())) ;
44
45 // Make pdf
46 RooRealVar tau("tau","tau",-1.54,-10,-0.1) ;
47 RooExponential model("model","model",t,tau) ;
48
49
50
51 // C r e a t e i n p u t d a t a
52 // ------------------------------------
53
54 // Generate complete dataset without acceptance cuts (for reference)
55 RooDataSet* dall = model.generate(t,10000) ;
56
57 // Generate a (fake) prototype dataset for acceptance limit values
58 RooDataSet* tmp = RooGaussian("gmin","gmin",tmin,RooConst(0),RooConst(0.5)).generate(tmin,5000) ;
59
60 // Generate dataset with t values that observe (t>tmin)
61 RooDataSet* dacc = model.generate(t,ProtoData(*tmp)) ;
62
63
64
65 // 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
66 // -----------------------------------------------------------------------
67
68 RooFitResult* r = model.fitTo(*dacc,Save()) ;
69
70
71
72 // 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
73 // ---------------------------------------------------------------------------------
74
75 // Make plot frame, add datasets and overlay model
76 RooPlot* frame = t.frame(Title("Fit to data with per-event acceptance")) ;
77 dall->plotOn(frame,MarkerColor(kRed),LineColor(kRed)) ;
78 model.plotOn(frame) ;
79 dacc->plotOn(frame) ;
80
81 // Print fit results to demonstrate absence of bias
82 r->Print("v") ;
83
84
85 new TCanvas("rf314_paramranges","rf314_paramranges",600,600) ;
86 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
87
88 return ;
89}
ROOT::R::TRInterface & r
Definition: Object.C:4
@ kRed
Definition: Rtypes.h:63
#define gPad
Definition: TVirtualPad.h:286
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:531
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition: RooAbsPdf.h:56
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Exponential p.d.f.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:550
RooCmdArg ProtoData(const RooDataSet &protoData, Bool_t randomizeOrder=kFALSE, Bool_t resample=kFALSE)
RooCmdArg MarkerColor(Color_t color)
RooConstVar & RooConst(Double_t val)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg LineColor(Color_t color)
const char * Title
Definition: TXMLSetup.cxx:67