Logo ROOT   6.14/05
Reference Guide
rf201_composite.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #201
5 ///
6 /// Composite p.d.f with signal and background component
7 ///
8 /// pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
9 ///
10 /// \macro_image
11 /// \macro_output
12 /// \macro_code
13 /// \author 07/2008 - Wouter Verkerke
14 
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "RooChebychev.h"
19 #include "RooAddPdf.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 using namespace RooFit ;
24 
25 
26 void rf201_composite()
27 {
28  // S e t u p c o m p o n e n t p d f s
29  // ---------------------------------------
30 
31  // Declare observable x
32  RooRealVar x("x","x",0,10) ;
33 
34  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
35  RooRealVar mean("mean","mean of gaussians",5) ;
36  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
37  RooRealVar sigma2("sigma2","width of gaussians",1) ;
38 
39  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
40  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
41 
42  // Build Chebychev polynomial p.d.f.
43  RooRealVar a0("a0","a0",0.5,0.,1.) ;
44  RooRealVar a1("a1","a1",0.2,0.,1.) ;
45  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
46 
47 
48  // ---------------------------------------------
49  // M E T H O D 1 - T w o R o o A d d P d f s
50  // =============================================
51 
52 
53  // A d d s i g n a l c o m p o n e n t s
54  // ------------------------------------------
55 
56  // Sum the signal components into a composite signal p.d.f.
57  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
58  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
59 
60 
61  // A d d s i g n a l a n d b a c k g r o u n d
62  // ------------------------------------------------
63 
64  // Sum the composite signal and background
65  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
66  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
67 
68 
69  // S a m p l e , f i t a n d p l o t m o d e l
70  // ---------------------------------------------------
71 
72  // Generate a data sample of 1000 events in x from model
73  RooDataSet *data = model.generate(x,1000) ;
74 
75  // Fit model to data
76  model.fitTo(*data) ;
77 
78  // Plot data and PDF overlaid
79  RooPlot* xframe = x.frame(Title("Example of composite pdf=(sig1+sig2)+bkg")) ;
80  data->plotOn(xframe) ;
81  model.plotOn(xframe) ;
82 
83  // Overlay the background component of model with a dashed line
84  model.plotOn(xframe,Components(bkg),LineStyle(kDashed)) ;
85 
86  // Overlay the background+sig2 components of model with a dotted line
87  model.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineStyle(kDotted)) ;
88 
89  // Print structure of composite p.d.f.
90  model.Print("t") ;
91 
92 
93  // ---------------------------------------------------------------------------------------------
94  // M E T H O D 2 - O n e R o o A d d P d f w i t h r e c u r s i v e f r a c t i o n s
95  // =============================================================================================
96 
97  // Construct sum of models on one go using recursive fraction interpretations
98  //
99  // model2 = bkg + (sig1 + sig2)
100  //
101  RooAddPdf model2("model","g1+g2+a",RooArgList(bkg,sig1,sig2),RooArgList(bkgfrac,sig1frac),kTRUE) ;
102 
103  // NB: Each coefficient is interpreted as the fraction of the
104  // left-hand component of the i-th recursive sum, i.e.
105  //
106  // sum4 = A + ( B + ( C + D) with fraction fA, fB and fC expands to
107  //
108  // sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
109 
110 
111  // P l o t r e c u r s i v e a d d i t i o n m o d e l
112  // ---------------------------------------------------------
113  model2.plotOn(xframe,LineColor(kRed),LineStyle(kDashed)) ;
114  model2.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineColor(kRed),LineStyle(kDashed)) ;
115  model2.Print("t") ;
116 
117 
118  // Draw the frame on the canvas
119  new TCanvas("rf201_composite","rf201_composite",600,600) ;
120  gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;
121 
122 
123 }
124 
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
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
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
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
RooCmdArg LineStyle(Style_t style)
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
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
RooCmdArg Components(const RooArgSet &compSet)
#define gPad
Definition: TVirtualPad.h:285
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558