Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf306_condpereventerrors.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// Multidimensional models: conditional p.d.f. with per-event errors
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooGaussModel.h"
16 #include "RooDecay.h"
17 #include "RooLandau.h"
18 #include "RooPlot.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "TH2D.h"
22 using namespace RooFit;
23 
24 void rf306_condpereventerrors()
25 {
26  // 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
27  // ----------------------------------------------------------------------------------------------
28 
29  // Observables
30  RooRealVar dt("dt", "dt", -10, 10);
31  RooRealVar dterr("dterr", "per-event error on dt", 0.01, 10);
32 
33  // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
34  RooRealVar bias("bias", "bias", 0, -10, 10);
35  RooRealVar sigma("sigma", "per-event error scale factor", 1, 0.1, 10);
36  RooGaussModel gm("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr);
37 
38  // Construct decay(dt) (x) gauss1(dt|dterr)
39  RooRealVar tau("tau", "tau", 1.548);
40  RooDecay decay_gm("decay_gm", "decay", dt, tau, gm, RooDecay::DoubleSided);
41 
42  // C o n s t r u c t f a k e ' e x t e r n a l ' d a t a w i t h p e r - e v e n t e r r o r
43  // ------------------------------------------------------------------------------------------------------
44 
45  // Use landau p.d.f to get somewhat realistic distribution with long tail
46  RooLandau pdfDtErr("pdfDtErr", "pdfDtErr", dterr, RooConst(1), RooConst(0.25));
47  RooDataSet *expDataDterr = pdfDtErr.generate(dterr, 10000);
48 
49  // S a m p l e d a t a f r o m c o n d i t i o n a l d e c a y _ g m ( d t | d t e r r )
50  // ---------------------------------------------------------------------------------------------
51 
52  // Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
53  RooDataSet *data = decay_gm.generate(dt, ProtoData(*expDataDterr));
54 
55  // 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 )
56  // ---------------------------------------------------------------------
57 
58  // Specify dterr as conditional observable
59  decay_gm.fitTo(*data, ConditionalObservables(dterr));
60 
61  // 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 )
62  // ---------------------------------------------------------------------
63 
64  // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
65  TH1 *hh_decay = decay_gm.createHistogram("hh_decay", dt, Binning(50), YVar(dterr, Binning(50)));
66  hh_decay->SetLineColor(kBlue);
67 
68  // Plot decay_gm(dt|dterr) at various values of dterr
69  RooPlot *frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr"));
70  for (Int_t ibin = 0; ibin < 100; ibin += 20) {
71  dterr.setBin(ibin);
72  decay_gm.plotOn(frame, Normalization(5.));
73  }
74 
75  // Make projection of data an dt
76  RooPlot *frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt"));
77  data->plotOn(frame2);
78 
79  // Make projection of decay(dt|dterr) on dt.
80  //
81  // Instead of integrating out dterr, make a weighted average of curves
82  // at values dterr_i as given in the external dataset.
83  // (The kTRUE argument bins the data before projection to speed up the process)
84  decay_gm.plotOn(frame2, ProjWData(*expDataDterr, kTRUE));
85 
86  // Draw all frames on canvas
87  TCanvas *c = new TCanvas("rf306_condpereventerrors", "rf306_condperventerrors", 1200, 400);
88  c->Divide(3);
89  c->cd(1);
90  gPad->SetLeftMargin(0.20);
91  hh_decay->GetZaxis()->SetTitleOffset(2.5);
92  hh_decay->Draw("surf");
93  c->cd(2);
94  gPad->SetLeftMargin(0.15);
95  frame->GetYaxis()->SetTitleOffset(1.6);
96  frame->Draw();
97  c->cd(3);
98  gPad->SetLeftMargin(0.15);
99  frame2->GetYaxis()->SetTitleOffset(1.6);
100  frame2->Draw();
101 }