Logo ROOT   6.14/05
Reference Guide
rf704_amplitudefit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook
4 /// 'SPECIAL PDFS' RooFit tutorial macro #704
5 ///
6 /// Using a p.d.f defined by a sum of real-valued amplitude components
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "RooTruthModel.h"
19 #include "RooFormulaVar.h"
20 #include "RooRealSumPdf.h"
21 #include "RooPolyVar.h"
22 #include "RooProduct.h"
23 #include "TH1.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "RooPlot.h"
27 using namespace RooFit ;
28 
29 
30 void rf704_amplitudefit()
31 {
32  // 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
33  // -------------------------------------------------------
34 
35  // Observables
36  RooRealVar t("t","time",-1.,15.);
37  RooRealVar cosa("cosa","cos(alpha)",-1.,1.);
38 
39  // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
40  RooRealVar tau("tau","#tau",1.5);
41  RooRealVar deltaGamma("deltaGamma","deltaGamma", 0.3);
42  RooTruthModel truthModel("tm","tm",t) ;
43  RooFormulaVar coshGBasis("coshGBasis","exp(-@0/ @1)*cosh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
44  RooFormulaVar sinhGBasis("sinhGBasis","exp(-@0/ @1)*sinh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
45  RooAbsReal* coshGConv = truthModel.convolution(&coshGBasis,&t);
46  RooAbsReal* sinhGConv = truthModel.convolution(&sinhGBasis,&t);
47 
48  // Construct polynomial amplitudes in cos(a)
49  RooPolyVar poly1("poly1","poly1",cosa,RooArgList(RooConst(0.5),RooConst(0.2),RooConst(0.2)),0);
50  RooPolyVar poly2("poly2","poly2",cosa,RooArgList(RooConst(1),RooConst(-0.2),RooConst(3)),0);
51 
52  // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
53  RooProduct ampl1("ampl1","amplitude 1",RooArgSet(poly1,*coshGConv));
54  RooProduct ampl2("ampl2","amplitude 2",RooArgSet(poly2,*sinhGConv));
55 
56 
57 
58  // C o n s t r u c t a m p l i t u d e s u m p d f
59  // -----------------------------------------------------
60 
61  // Amplitude strengths
62  RooRealVar f1("f1","f1",1,0,2) ;
63  RooRealVar f2("f2","f2",0.5,0,2) ;
64 
65  // Construct pdf
66  RooRealSumPdf pdf("pdf","pdf",RooArgList(ampl1,ampl2),RooArgList(f1,f2)) ;
67 
68  // Generate some toy data from pdf
69  RooDataSet* data = pdf.generate(RooArgSet(t,cosa),10000);
70 
71  // Fit pdf to toy data with only amplitude strength floating
72  pdf.fitTo(*data) ;
73 
74 
75 
76  // P l o t a m p l i t u d e s u m p d f
77  // -------------------------------------------
78 
79  // Make 2D plots of amplitudes
80  TH1* hh_cos = ampl1.createHistogram("hh_cos",t,Binning(50),YVar(cosa,Binning(50))) ;
81  TH1* hh_sin = ampl2.createHistogram("hh_sin",t,Binning(50),YVar(cosa,Binning(50))) ;
82  hh_cos->SetLineColor(kBlue) ;
83  hh_sin->SetLineColor(kRed) ;
84 
85 
86  // Make projection on t, plot data, pdf and its components
87  // Note component projections may be larger than sum because amplitudes can be negative
88  RooPlot* frame1 = t.frame();
89  data->plotOn(frame1);
90  pdf.plotOn(frame1);
91  pdf.plotOn(frame1,Components(ampl1),LineStyle(kDashed));
92  pdf.plotOn(frame1,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
93 
94  // Make projection on cosa, plot data, pdf and its components
95  // Note that components projection may be larger than sum because amplitudes can be negative
96  RooPlot* frame2 = cosa.frame();
97  data->plotOn(frame2);
98  pdf.plotOn(frame2);
99  pdf.plotOn(frame2,Components(ampl1),LineStyle(kDashed));
100  pdf.plotOn(frame2,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
101 
102 
103 
104  TCanvas* c = new TCanvas("rf704_amplitudefit","rf704_amplitudefit",800,800) ;
105  c->Divide(2,2) ;
106  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw();
107  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw();
108  c->cd(3) ; gPad->SetLeftMargin(0.20) ; hh_cos->GetZaxis()->SetTitleOffset(2.3) ; hh_cos->Draw("surf") ;
109  c->cd(4) ; gPad->SetLeftMargin(0.20) ; hh_sin->GetZaxis()->SetTitleOffset(2.3) ; hh_sin->Draw("surf") ;
110 
111 
112 }
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1117
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:568
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:59
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
Class RooRealSumPdf implements a PDF constructed from a sum of functions:
Definition: RooRealSumPdf.h:24
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
RooCmdArg LineStyle(Style_t style)
Class RooPolyVar is a RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficie...
Definition: RooPolyVar.h:28
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooProduct a RooAbsReal implementation that represent the product of a given set of other RooAbsReal ...
Definition: RooProduct.h:32
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The Canvas class.
Definition: TCanvas.h:31
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetZaxis()
Definition: TH1.h:317
RooConstVar & RooConst(Double_t val)
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
TF1 * f1
Definition: legend1.C:11
#define gPad
Definition: TVirtualPad.h:285
#define c(i)
Definition: RSha256.hxx:101
Definition: Rtypes.h:59
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooCmdArg Binning(const RooAbsBinning &binning)