Logo ROOT  
Reference Guide
rf314_paramfitrange.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4///
5/// Multidimensional models: working with parametrized ranges in a fit.
6/// This an example of a fit with an acceptance that changes per-event
7///
8/// `pdf = exp(-t/tau)` with `t[tmin,5]`
9///
10/// where `t` and `tmin` are both observables in the dataset
11///
12/// \macro_image
13/// \macro_output
14/// \macro_code
15///
16/// \date 07/2008
17/// \author Wouter Verkerke
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
32{
33
34 // 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
35 // ---------------------------------------------------------------
36
37 // Declare observables
38 RooRealVar t("t", "t", 0, 5);
39 RooRealVar tmin("tmin", "tmin", 0, 0, 5);
40
41 // Make parametrized range in t : [tmin,5]
42 t.setRange(tmin, RooConst(t.getMax()));
43
44 // Make pdf
45 RooRealVar tau("tau", "tau", -1.54, -10, -0.1);
46 RooExponential model("model", "model", t, tau);
47
48 // C r e a t e i n p u t d a t a
49 // ------------------------------------
50
51 // Generate complete dataset without acceptance cuts (for reference)
52 RooDataSet *dall = model.generate(t, 10000);
53
54 // Generate a (fake) prototype dataset for acceptance limit values
55 RooDataSet *tmp = RooGaussian("gmin", "gmin", tmin, RooConst(0), RooConst(0.5)).generate(tmin, 5000);
56
57 // Generate dataset with t values that observe (t>tmin)
58 RooDataSet *dacc = model.generate(t, ProtoData(*tmp));
59
60 // 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
61 // -----------------------------------------------------------------------
62
63 RooFitResult *r = model.fitTo(*dacc, Save());
64
65 // 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
66 // ---------------------------------------------------------------------------------
67
68 // Make plot frame, add datasets and overlay model
69 RooPlot *frame = t.frame(Title("Fit to data with per-event acceptance"));
70 dall->plotOn(frame, MarkerColor(kRed), LineColor(kRed));
71 model.plotOn(frame);
72 dacc->plotOn(frame);
73
74 // Print fit results to demonstrate absence of bias
75 r->Print("v");
76
77 new TCanvas("rf314_paramranges", "rf314_paramranges", 600, 600);
78 gPad->SetLeftMargin(0.15);
79 frame->GetYaxis()->SetTitleOffset(1.6);
80 frame->Draw();
81
82 return;
83}
ROOT::R::TRInterface & r
Definition: Object.C:4
@ kRed
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:287
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:546
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:55
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
Exponential PDF.
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:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition: RooPlot.cxx:277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:27
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:550
RooConstVar & RooConst(Double_t val)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg ProtoData(const RooDataSet &protoData, Bool_t randomizeOrder=kFALSE, Bool_t resample=kFALSE)
RooCmdArg MarkerColor(Color_t color)
RooCmdArg LineColor(Color_t color)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition: TXMLSetup.cxx:67