Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf704_amplitudefit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// Special pdf's: using a pdf defined by a sum of real-valued amplitude components
5///
6/// \macro_image
7/// \macro_output
8/// \macro_code
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooConstVar.h"
17#include "RooTruthModel.h"
18#include "RooFormulaVar.h"
19#include "RooRealSumPdf.h"
20#include "RooPolyVar.h"
21#include "RooProduct.h"
22#include "TH1.h"
23#include "TCanvas.h"
24#include "TAxis.h"
25#include "RooPlot.h"
26using namespace RooFit;
27
29{
30 // 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
31 // -------------------------------------------------------
32
33 // Observables
34 RooRealVar t("t", "time", -1., 15.);
35 RooRealVar cosa("cosa", "cos(alpha)", -1., 1.);
36
37 // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
38 RooRealVar tau("tau", "#tau", 1.5);
39 RooRealVar deltaGamma("deltaGamma", "deltaGamma", 0.3);
40 RooTruthModel truthModel("tm", "tm", t);
41 RooFormulaVar coshGBasis("coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
42 RooFormulaVar sinhGBasis("sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
43 RooAbsReal *coshGConv = truthModel.convolution(&coshGBasis, &t);
44 RooAbsReal *sinhGConv = truthModel.convolution(&sinhGBasis, &t);
45
46 // Construct polynomial amplitudes in cos(a)
47 RooPolyVar poly1("poly1", "poly1", cosa, RooArgList(RooConst(0.5), RooConst(0.2), RooConst(0.2)), 0);
48 RooPolyVar poly2("poly2", "poly2", cosa, RooArgList(RooConst(1), RooConst(-0.2), RooConst(3)), 0);
49
50 // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
51 RooProduct ampl1("ampl1", "amplitude 1", RooArgSet(poly1, *coshGConv));
52 RooProduct ampl2("ampl2", "amplitude 2", RooArgSet(poly2, *sinhGConv));
53
54 // C o n s t r u c t a m p l i t u d e s u m p d f
55 // -----------------------------------------------------
56
57 // Amplitude strengths
58 RooRealVar f1("f1", "f1", 1, 0, 2);
59 RooRealVar f2("f2", "f2", 0.5, 0, 2);
60
61 // Construct pdf
62 RooRealSumPdf pdf("pdf", "pdf", RooArgList(ampl1, ampl2), RooArgList(f1, f2));
63
64 // Generate some toy data from pdf
65 RooDataSet *data = pdf.generate(RooArgSet(t, cosa), 10000);
66
67 // Fit pdf to toy data with only amplitude strength floating
68 pdf.fitTo(*data);
69
70 // P l o t a m p l i t u d e s u m p d f
71 // -------------------------------------------
72
73 // Make 2D plots of amplitudes
74 TH1 *hh_cos = ampl1.createHistogram("hh_cos", t, Binning(50), YVar(cosa, Binning(50)));
75 TH1 *hh_sin = ampl2.createHistogram("hh_sin", t, Binning(50), YVar(cosa, Binning(50)));
76 hh_cos->SetLineColor(kBlue);
77 hh_sin->SetLineColor(kRed);
78
79 // Make projection on t, plot data, pdf and its components
80 // Note component projections may be larger than sum because amplitudes can be negative
81 RooPlot *frame1 = t.frame();
82 data->plotOn(frame1);
83 pdf.plotOn(frame1);
84 pdf.plotOn(frame1, Components(ampl1), LineStyle(kDashed));
85 pdf.plotOn(frame1, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
86
87 // Make projection on cosa, plot data, pdf and its components
88 // Note that components projection may be larger than sum because amplitudes can be negative
89 RooPlot *frame2 = cosa.frame();
90 data->plotOn(frame2);
91 pdf.plotOn(frame2);
92 pdf.plotOn(frame2, Components(ampl1), LineStyle(kDashed));
93 pdf.plotOn(frame2, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
94
95 TCanvas *c = new TCanvas("rf704_amplitudefit", "rf704_amplitudefit", 800, 800);
96 c->Divide(2, 2);
97 c->cd(1);
98 gPad->SetLeftMargin(0.15);
99 frame1->GetYaxis()->SetTitleOffset(1.4);
100 frame1->Draw();
101 c->cd(2);
102 gPad->SetLeftMargin(0.15);
103 frame2->GetYaxis()->SetTitleOffset(1.4);
104 frame2->Draw();
105 c->cd(3);
106 gPad->SetLeftMargin(0.20);
107 hh_cos->GetZaxis()->SetTitleOffset(2.3);
108 hh_cos->Draw("surf");
109 c->cd(4);
110 gPad->SetLeftMargin(0.20);
111 hh_sin->GetZaxis()->SetTitleOffset(2.3);
112 hh_sin->Draw("surf");
113}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
#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
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
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
Class RooPolyVar is a RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficie...
Definition RooPolyVar.h:28
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
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 Binning(const RooAbsBinning &binning)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
TF1 * f1
Definition legend1.C:11
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...