Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf307_fullpereventerrors.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// Multidimensional models: full p.d.f. with per-event errors
5 ///
6 ///
7 ///
8 /// \macro_code
9 /// \author 07/2008 - 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"
24 using namespace RooFit;
25 
26 void rf307_fullpereventerrors()
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 p.d.f 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 p.d.f.
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 p.d.f 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 }