ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf201_composite.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #201
4 //
5 // Composite p.d.f with signal and background component
6 //
7 // pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
8 //
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooDataSet.h"
19 #include "RooGaussian.h"
20 #include "RooChebychev.h"
21 #include "RooAddPdf.h"
22 #include "TCanvas.h"
23 #include "RooPlot.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic201 : public RooFitTestUnit
28 {
29 public:
30  TestBasic201(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Addition operator p.d.f.",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // S e t u p c o m p o n e n t p d f s
34  // ---------------------------------------
35 
36  // Declare observable x
37  RooRealVar x("x","x",0,10) ;
38 
39  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
40  RooRealVar mean("mean","mean of gaussians",5) ;
41  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
42  RooRealVar sigma2("sigma2","width of gaussians",1) ;
43 
44  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
45  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
46 
47  // Build Chebychev polynomial p.d.f.
48  RooRealVar a0("a0","a0",0.5,0.,1.) ;
49  RooRealVar a1("a1","a1",-0.2,0.,1.) ;
50  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
51 
52 
53  ////////////////////////////////////////////////////
54  // M E T H O D 1 - T w o R o o A d d P d f s //
55  ////////////////////////////////////////////////////
56 
57 
58  // A d d s i g n a l c o m p o n e n t s
59  // ------------------------------------------
60 
61  // Sum the signal components into a composite signal p.d.f.
62  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
63  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
64 
65 
66  // A d d s i g n a l a n d b a c k g r o u n d
67  // ------------------------------------------------
68 
69  // Sum the composite signal and background
70  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
71  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
72 
73 
74  // S a m p l e , f i t a n d p l o t m o d e l
75  // ---------------------------------------------------
76 
77  // Generate a data sample of 1000 events in x from model
78  RooDataSet *data = model.generate(x,1000) ;
79 
80  // Fit model to data
81  model.fitTo(*data) ;
82 
83  // Plot data and PDF overlaid
84  RooPlot* xframe = x.frame(Title("Example of composite pdf=(sig1+sig2)+bkg")) ;
85  data->plotOn(xframe) ;
86  model.plotOn(xframe) ;
87 
88  // Overlay the background component of model with a dashed line
89  model.plotOn(xframe,Components(bkg),LineStyle(kDashed)) ;
90 
91  // Overlay the background+sig2 components of model with a dotted line
92  model.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineStyle(kDotted)) ;
93 
94 
95  ////////////////////////////////////////////////////////////////////////////////////////////////////
96  // 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 //
97  ////////////////////////////////////////////////////////////////////////////////////////////////////
98 
99  // Construct sum of models on one go using recursive fraction interpretations
100  //
101  // model2 = bkg + (sig1 + sig2)
102  //
103  RooAddPdf model2("model","g1+g2+a",RooArgList(bkg,sig1,sig2),RooArgList(bkgfrac,sig1frac),kTRUE) ;
104 
105  // NB: Each coefficient is interpreted as the fraction of the
106  // left-hand component of the i-th recursive sum, i.e.
107  //
108  // sum4 = A + ( B + ( C + D) with fraction fA, fB and fC expands to
109  //
110  // sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
111 
112 
113  // 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
114  // ---------------------------------------------------------
115  model2.plotOn(xframe,LineColor(kRed),LineStyle(kDashed)) ;
116  model2.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineColor(kRed),LineStyle(kDashed)) ;
117 
118 
119  regPlot(xframe,"rf201_plot1") ;
120 
121  delete data ;
122  return kTRUE ;
123 
124  }
125 
126 } ;
127 
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:61
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.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
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:626
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
bool verbose
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
RooCmdArg Components(const RooArgSet &compSet)
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const Bool_t kTRUE
Definition: Rtypes.h:91