Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf703_effpdfprod.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Speecial p.d.f.'s: using a product of an (acceptance) efficiency and a p.d.f as p.d.f.
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 "RooExponential.h"
16 #include "RooEffProd.h"
17 #include "RooFormulaVar.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 using namespace RooFit;
22 
23 void rf703_effpdfprod()
24 {
25  // 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
26  // ---------------------------------------------------------------
27 
28  // Declare observables
29  RooRealVar t("t", "t", 0, 5);
30 
31  // Make pdf
32  RooRealVar tau("tau", "tau", -1.54, -4, -0.1);
33  RooExponential model("model", "model", t, tau);
34 
35  // D e f i n e e f f i c i e n c y f u n c t i o n
36  // ---------------------------------------------------
37 
38  // Use error function to simulate turn-on slope
39  RooFormulaVar eff("eff", "0.5*(TMath::Erf((t-1)/0.5)+1)", t);
40 
41  // D e f i n e d e c a y p d f w i t h e f f i c i e n c y
42  // ---------------------------------------------------------------
43 
44  // Multiply pdf(t) with efficiency in t
45  RooEffProd modelEff("modelEff", "model with efficiency", model, eff);
46 
47  // P l o t e f f i c i e n c y , p d f
48  // ----------------------------------------
49 
50  RooPlot *frame1 = t.frame(Title("Efficiency"));
51  eff.plotOn(frame1, LineColor(kRed));
52 
53  RooPlot *frame2 = t.frame(Title("Pdf with and without efficiency"));
54 
55  model.plotOn(frame2, LineStyle(kDashed));
56  modelEff.plotOn(frame2);
57 
58  // G e n e r a t e t o y d a t a , f i t m o d e l E f f t o d a t a
59  // ------------------------------------------------------------------------------
60 
61  // Generate events. If the input pdf has an internal generator, the internal generator
62  // is used and an accept/reject sampling on the efficiency is applied.
63  RooDataSet *data = modelEff.generate(t, 10000);
64 
65  // Fit pdf. The normalization integral is calculated numerically.
66  modelEff.fitTo(*data);
67 
68  // Plot generated data and overlay fitted pdf
69  RooPlot *frame3 = t.frame(Title("Fitted pdf with efficiency"));
70  data->plotOn(frame3);
71  modelEff.plotOn(frame3);
72 
73  TCanvas *c = new TCanvas("rf703_effpdfprod", "rf703_effpdfprod", 1200, 400);
74  c->Divide(3);
75  c->cd(1);
76  gPad->SetLeftMargin(0.15);
77  frame1->GetYaxis()->SetTitleOffset(1.4);
78  frame1->Draw();
79  c->cd(2);
80  gPad->SetLeftMargin(0.15);
81  frame2->GetYaxis()->SetTitleOffset(1.6);
82  frame2->Draw();
83  c->cd(3);
84  gPad->SetLeftMargin(0.15);
85  frame3->GetYaxis()->SetTitleOffset(1.6);
86  frame3->Draw();
87 }