Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf307_fullpereventerrors.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// Multidimensional models: full pdf with per-event errors
5///
6/// \macro_code
7///
8/// \date July 2008
9/// \author Wouter Verkerke
10
11#include "RooRealVar.h"
12#include "RooDataSet.h"
13#include "RooGaussian.h"
14#include "RooGaussModel.h"
15#include "RooConstVar.h"
16#include "RooDecay.h"
17#include "RooLandau.h"
18#include "RooProdPdf.h"
19#include "RooHistPdf.h"
20#include "RooPlot.h"
21#include "TCanvas.h"
22#include "TAxis.h"
23#include "TH1.h"
24using namespace RooFit;
25
27{
28 // 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
29 // ----------------------------------------------------------------------------------------------
30
31 // Observables
32 RooRealVar dt("dt", "dt", -10, 10);
33 RooRealVar dterr("dterr", "per-event error on dt", 0.01, 10);
34
35 // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
36 RooRealVar bias("bias", "bias", 0, -10, 10);
37 RooRealVar sigma("sigma", "per-event error scale factor", 1, 0.1, 10);
38 RooGaussModel gm("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr);
39
40 // Construct decay(dt) (x) gauss1(dt|dterr)
41 RooRealVar tau("tau", "tau", 1.548);
42 RooDecay decay_gm("decay_gm", "decay", dt, tau, gm, RooDecay::DoubleSided);
43
44 // 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
45 // -----------------------------------------------------------------
46
47 // Use landau pdf to get empirical distribution with long tail
48 RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, RooConst(1), RooConst(0.25));
49 RooDataSet *expDataDterr = pdfDtErr.generate(dterr, 10000);
50
51 // Construct a histogram pdf to describe the shape of the dtErr distribution
52 RooDataHist *expHistDterr = expDataDterr->binnedClone();
53 RooHistPdf pdfErr("pdfErr", "pdfErr", dterr, *expHistDterr);
54
55 // 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
56 // r )
57 // ----------------------------------------------------------------------------------------------------------------------
58
59 // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
60 RooProdPdf model("model", "model", pdfErr, Conditional(decay_gm, dt));
61
62 // (Alternatively you could also use the landau shape pdfDtErr)
63 // RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
64
65 // 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
66 // ------------------------------------------------------------------
67
68 // Specify external dataset with dterr values to use model_dm as conditional pdf
69 RooDataSet *data = model.generate(RooArgSet(dt, dterr), 10000);
70
71 // 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 )
72 // ---------------------------------------------------------------------
73
74 // Specify dterr as conditional observable
75 model.fitTo(*data);
76
77 // 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 )
78 // ---------------------------------------------------------------------
79
80 // Make two-dimensional plot of conditional pdf in (dt,dterr)
81 TH1 *hh_model = model.createHistogram("hh_model", dt, Binning(50), YVar(dterr, Binning(50)));
82 hh_model->SetLineColor(kBlue);
83
84 // Make projection of data an dt
85 RooPlot *frame = dt.frame(Title("Projection of model(dt|dterr) on dt"));
86 data->plotOn(frame);
87 model.plotOn(frame);
88
89 // Draw all frames on canvas
90 TCanvas *c = new TCanvas("rf307_fullpereventerrors", "rf307_fullperventerrors", 800, 400);
91 c->Divide(2);
92 c->cd(1);
93 gPad->SetLeftMargin(0.20);
94 hh_model->GetZaxis()->SetTitleOffset(2.5);
95 hh_model->Draw("surf");
96 c->cd(2);
97 gPad->SetLeftMargin(0.15);
98 frame->GetYaxis()->SetTitleOffset(1.6);
99 frame->Draw();
100}
#define c(i)
Definition RSha256.hxx:101
@ kBlue
Definition Rtypes.h:66
#define gPad
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
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:37
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.
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.
RooHistPdf implements a probablity density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:29
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:1263
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:249
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:37
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:293
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:322
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
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:68