Logo ROOT   6.10/09
Reference Guide
rf704_amplitudefit.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'SPECIAL PDFS' RooFit tutorial macro #704
4 //
5 // Using a p.d.f defined by a sum of real-valued amplitude components
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooTruthModel.h"
20 #include "RooFormulaVar.h"
21 #include "RooRealSumPdf.h"
22 #include "RooPolyVar.h"
23 #include "RooProduct.h"
24 #include "TH1.h"
25 #include "TCanvas.h"
26 #include "RooPlot.h"
27 using namespace RooFit ;
28 
29 
30 // Elementary operations on a gaussian PDF
31 class TestBasic704 : public RooFitTestUnit
32 {
33 public:
34  TestBasic704(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Amplitude sum operator p.d.f",refFile,writeRef,verbose) {} ;
35  Bool_t testCode() {
36 
37  // 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
38  // -------------------------------------------------------
39 
40  // Observables
41  RooRealVar t("t","time",-1.,15.);
42  RooRealVar cosa("cosa","cos(alpha)",-1.,1.);
43 
44  // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
45  RooRealVar tau("tau","#tau",1.5);
46  RooRealVar deltaGamma("deltaGamma","deltaGamma", 0.3);
47  RooTruthModel tm("tm","tm",t) ;
48  RooFormulaVar coshGBasis("coshGBasis","exp(-@0/ @1)*cosh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
49  RooFormulaVar sinhGBasis("sinhGBasis","exp(-@0/ @1)*sinh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
50  RooAbsReal* coshGConv = tm.convolution(&coshGBasis,&t);
51  RooAbsReal* sinhGConv = tm.convolution(&sinhGBasis,&t);
52 
53  // Construct polynomial amplitudes in cos(a)
54  RooPolyVar poly1("poly1","poly1",cosa,RooArgList(RooConst(0.5),RooConst(0.2),RooConst(0.2)),0);
55  RooPolyVar poly2("poly2","poly2",cosa,RooArgList(RooConst(1),RooConst(-0.2),RooConst(3)),0);
56 
57  // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
58  RooProduct ampl1("ampl1","amplitude 1",RooArgSet(poly1,*coshGConv));
59  RooProduct ampl2("ampl2","amplitude 2",RooArgSet(poly2,*sinhGConv));
60 
61 
62 
63  // C o n s t r u c t a m p l i t u d e s u m p d f
64  // -----------------------------------------------------
65 
66  // Amplitude strengths
67  RooRealVar f1("f1","f1",1,0,2) ;
68  RooRealVar f2("f2","f2",0.5,0,2) ;
69 
70  // Construct pdf
71  RooRealSumPdf pdf("pdf","pdf",RooArgList(ampl1,ampl2),RooArgList(f1,f2)) ;
72 
73  // Generate some toy data from pdf
74  RooDataSet* data = pdf.generate(RooArgSet(t,cosa),10000);
75 
76  // Fit pdf to toy data with only amplitude strength floating
77  pdf.fitTo(*data) ;
78 
79 
80 
81  // P l o t a m p l i t u d e s u m p d f
82  // -------------------------------------------
83 
84  // Make 2D plots of amplitudes
85  TH1* hh_cos = ampl1.createHistogram("hh_cos",t,Binning(50),YVar(cosa,Binning(50))) ;
86  TH1* hh_sin = ampl2.createHistogram("hh_sin",t,Binning(50),YVar(cosa,Binning(50))) ;
87  hh_cos->SetLineColor(kBlue) ;
88  hh_sin->SetLineColor(kBlue) ;
89 
90 
91  // Make projection on t, plot data, pdf and its components
92  // Note component projections may be larger than sum because amplitudes can be negative
93  RooPlot* frame1 = t.frame();
94  data->plotOn(frame1);
95  pdf.plotOn(frame1);
96  pdf.plotOn(frame1,Components(ampl1),LineStyle(kDashed));
97  pdf.plotOn(frame1,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
98 
99  // Make projection on cosa, plot data, pdf and its components
100  // Note that components projection may be larger than sum because amplitudes can be negative
101  RooPlot* frame2 = cosa.frame();
102  data->plotOn(frame2);
103  pdf.plotOn(frame2);
104  pdf.plotOn(frame2,Components(ampl1),LineStyle(kDashed));
105  pdf.plotOn(frame2,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
106 
107 
108  regPlot(frame1,"rf704_plot1") ;
109  regPlot(frame2,"rf704_plot2") ;
110  regTH(hh_cos,"rf704_hh_cos") ;
111  regTH(hh_sin,"rf704_hh_sin") ;
112 
113  delete data ;
114  delete coshGConv ;
115  delete sinhGConv ;
116 
117  return kTRUE ;
118  }
119 } ;
double poly2(const double *x, const double *p)
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
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:552
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
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
bool verbose
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
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
RooConstVar & RooConst(Double_t val)
double f2(const double *x)
TF1 * f1
Definition: legend1.C:11
Definition: Rtypes.h:56
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooCmdArg Binning(const RooAbsBinning &binning)