Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf704_amplitudefit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// Speecial p.d.f.'s: using a p.d.f defined by a sum of real-valued amplitude components
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 "RooTruthModel.h"
16 #include "RooFormulaVar.h"
17 #include "RooRealSumPdf.h"
18 #include "RooPolyVar.h"
19 #include "RooProduct.h"
20 #include "TH1.h"
21 #include "TCanvas.h"
22 #include "TAxis.h"
23 #include "RooPlot.h"
24 using namespace RooFit;
25 
26 void rf704_amplitudefit()
27 {
28  // S e t u p 2 D a m p l i t u d e f u n c t i o n s
29  // -------------------------------------------------------
30 
31  // Observables
32  RooRealVar t("t", "time", -1., 15.);
33  RooRealVar cosa("cosa", "cos(alpha)", -1., 1.);
34 
35  // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
36  RooRealVar tau("tau", "#tau", 1.5);
37  RooRealVar deltaGamma("deltaGamma", "deltaGamma", 0.3);
38  RooTruthModel truthModel("tm", "tm", t);
39  RooFormulaVar coshGBasis("coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
40  RooFormulaVar sinhGBasis("sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
41  RooAbsReal *coshGConv = truthModel.convolution(&coshGBasis, &t);
42  RooAbsReal *sinhGConv = truthModel.convolution(&sinhGBasis, &t);
43 
44  // Construct polynomial amplitudes in cos(a)
45  RooPolyVar poly1("poly1", "poly1", cosa, RooArgList(RooConst(0.5), RooConst(0.2), RooConst(0.2)), 0);
46  RooPolyVar poly2("poly2", "poly2", cosa, RooArgList(RooConst(1), RooConst(-0.2), RooConst(3)), 0);
47 
48  // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
49  RooProduct ampl1("ampl1", "amplitude 1", RooArgSet(poly1, *coshGConv));
50  RooProduct ampl2("ampl2", "amplitude 2", RooArgSet(poly2, *sinhGConv));
51 
52  // C o n s t r u c t a m p l i t u d e s u m p d f
53  // -----------------------------------------------------
54 
55  // Amplitude strengths
56  RooRealVar f1("f1", "f1", 1, 0, 2);
57  RooRealVar f2("f2", "f2", 0.5, 0, 2);
58 
59  // Construct pdf
60  RooRealSumPdf pdf("pdf", "pdf", RooArgList(ampl1, ampl2), RooArgList(f1, f2));
61 
62  // Generate some toy data from pdf
63  RooDataSet *data = pdf.generate(RooArgSet(t, cosa), 10000);
64 
65  // Fit pdf to toy data with only amplitude strength floating
66  pdf.fitTo(*data);
67 
68  // P l o t a m p l i t u d e s u m p d f
69  // -------------------------------------------
70 
71  // Make 2D plots of amplitudes
72  TH1 *hh_cos = ampl1.createHistogram("hh_cos", t, Binning(50), YVar(cosa, Binning(50)));
73  TH1 *hh_sin = ampl2.createHistogram("hh_sin", t, Binning(50), YVar(cosa, Binning(50)));
74  hh_cos->SetLineColor(kBlue);
75  hh_sin->SetLineColor(kRed);
76 
77  // Make projection on t, plot data, pdf and its components
78  // Note component projections may be larger than sum because amplitudes can be negative
79  RooPlot *frame1 = t.frame();
80  data->plotOn(frame1);
81  pdf.plotOn(frame1);
82  pdf.plotOn(frame1, Components(ampl1), LineStyle(kDashed));
83  pdf.plotOn(frame1, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
84 
85  // Make projection on cosa, plot data, pdf and its components
86  // Note that components projection may be larger than sum because amplitudes can be negative
87  RooPlot *frame2 = cosa.frame();
88  data->plotOn(frame2);
89  pdf.plotOn(frame2);
90  pdf.plotOn(frame2, Components(ampl1), LineStyle(kDashed));
91  pdf.plotOn(frame2, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
92 
93  TCanvas *c = new TCanvas("rf704_amplitudefit", "rf704_amplitudefit", 800, 800);
94  c->Divide(2, 2);
95  c->cd(1);
96  gPad->SetLeftMargin(0.15);
97  frame1->GetYaxis()->SetTitleOffset(1.4);
98  frame1->Draw();
99  c->cd(2);
100  gPad->SetLeftMargin(0.15);
101  frame2->GetYaxis()->SetTitleOffset(1.4);
102  frame2->Draw();
103  c->cd(3);
104  gPad->SetLeftMargin(0.20);
105  hh_cos->GetZaxis()->SetTitleOffset(2.3);
106  hh_cos->Draw("surf");
107  c->cd(4);
108  gPad->SetLeftMargin(0.20);
109  hh_sin->GetZaxis()->SetTitleOffset(2.3);
110  hh_sin->Draw("surf");
111 }