Logo ROOT  
Reference Guide
rf307_fullpereventerrors.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4///
5/// Multidimensional models: full p.d.f. with per-event errors
6///
7/// \macro_code
8///
9/// \date 07/2008
10/// \author Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooDataSet.h"
14#include "RooGaussian.h"
15#include "RooGaussModel.h"
16#include "RooConstVar.h"
17#include "RooDecay.h"
18#include "RooLandau.h"
19#include "RooProdPdf.h"
20#include "RooHistPdf.h"
21#include "RooPlot.h"
22#include "TCanvas.h"
23#include "TAxis.h"
24#include "TH1.h"
25using namespace RooFit;
26
28{
29 // B - p h y s i c s p d f w i t h p e r - e v e n t G a u s s i a n r e s o l u t i o n
30 // ----------------------------------------------------------------------------------------------
31
32 // Observables
33 RooRealVar dt("dt", "dt", -10, 10);
34 RooRealVar dterr("dterr", "per-event error on dt", 0.01, 10);
35
36 // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
37 RooRealVar bias("bias", "bias", 0, -10, 10);
38 RooRealVar sigma("sigma", "per-event error scale factor", 1, 0.1, 10);
39 RooGaussModel gm("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr);
40
41 // Construct decay(dt) (x) gauss1(dt|dterr)
42 RooRealVar tau("tau", "tau", 1.548);
43 RooDecay decay_gm("decay_gm", "decay", dt, tau, gm, RooDecay::DoubleSided);
44
45 // C o n s t r u c t e m p i r i c a l p d f f o r p e r - e v e n t e r r o r
46 // -----------------------------------------------------------------
47
48 // Use landau p.d.f to get empirical distribution with long tail
49 RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, RooConst(1), RooConst(0.25));
50 RooDataSet *expDataDterr = pdfDtErr.generate(dterr, 10000);
51
52 // Construct a histogram pdf to describe the shape of the dtErr distribution
53 RooDataHist *expHistDterr = expDataDterr->binnedClone();
54 RooHistPdf pdfErr("pdfErr", "pdfErr", dterr, *expHistDterr);
55
56 // C o n s t r u c t c o n d i t i o n a l p r o d u c t d e c a y _ d m ( d t | d t e r r ) * p d f ( d t e r
57 // r )
58 // ----------------------------------------------------------------------------------------------------------------------
59
60 // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
61 RooProdPdf model("model", "model", pdfErr, Conditional(decay_gm, dt));
62
63 // (Alternatively you could also use the landau shape pdfDtErr)
64 // RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
65
66 // S a m p l e, f i t a n d p l o t p r o d u c t m o d e l
67 // ------------------------------------------------------------------
68
69 // Specify external dataset with dterr values to use model_dm as conditional p.d.f.
70 RooDataSet *data = model.generate(RooArgSet(dt, dterr), 10000);
71
72 // F i t c o n d i t i o n a l d e c a y _ d m ( d t | d t e r r )
73 // ---------------------------------------------------------------------
74
75 // Specify dterr as conditional observable
76 model.fitTo(*data);
77
78 // P l o t c o n d i t i o n a l d e c a y _ d m ( d t | d t e r r )
79 // ---------------------------------------------------------------------
80
81 // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
82 TH1 *hh_model = model.createHistogram("hh_model", dt, Binning(50), YVar(dterr, Binning(50)));
83 hh_model->SetLineColor(kBlue);
84
85 // Make projection of data an dt
86 RooPlot *frame = dt.frame(Title("Projection of model(dt|dterr) on dt"));
87 data->plotOn(frame);
88 model.plotOn(frame);
89
90 // Draw all frames on canvas
91 TCanvas *c = new TCanvas("rf307_fullpereventerrors", "rf307_fullperventerrors", 800, 400);
92 c->Divide(2);
93 c->cd(1);
94 gPad->SetLeftMargin(0.20);
95 hh_model->GetZaxis()->SetTitleOffset(2.5);
96 hh_model->Draw("surf");
97 c->cd(2);
98 gPad->SetLeftMargin(0.15);
99 frame->GetYaxis()->SetTitleOffset(1.6);
100 frame->Draw();
101}
#define c(i)
Definition: RSha256.hxx:101
@ kBlue
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
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
RooDataHist * binnedClone(const char *newName=0, const char *newTitle=0) const
Return binned clone of this dataset.
Definition: RooDataSet.cxx:951
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
@ DoubleSided
Definition: RooDecay.h:25
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
RooHistPdf implements a probablity density function sampled from a multidimensional histogram.
Definition: RooHistPdf.h:28
Landau distribution p.d.f.
Definition: RooLandau.h:24
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
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
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
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:27
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooConstVar & RooConst(Double_t val)
RooCmdArg Conditional(const RooArgSet &pdfSet, const RooArgSet &depSet, Bool_t depsAreCond=kFALSE)
RooCmdArg Binning(const RooAbsBinning &binning)
const Double_t sigma
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition: TXMLSetup.cxx:67