Logo ROOT   6.10/09
Reference Guide
rf207_comptools.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #207
5 ///
6 /// Tools and utilities for manipulation of composite objects
7 ///
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "RooChebychev.h"
17 #include "RooExponential.h"
18 #include "RooAddPdf.h"
19 #include "RooPlot.h"
20 #include "RooCustomizer.h"
21 #include "TCanvas.h"
22 #include "TAxis.h"
23 #include "TH1.h"
24 using namespace RooFit ;
25 
26 
27 
28 void rf207_comptools()
29 {
30 
31  // S e t u p c o m p o s i t e p d f, d a t a s e t
32  // --------------------------------------------------------
33 
34  // Declare observable x
35  RooRealVar x("x","x",0,10) ;
36 
37  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
38  RooRealVar mean("mean","mean of gaussians",5) ;
39  RooRealVar sigma("sigma","width of gaussians",0.5) ;
40  RooGaussian sig("sig","Signal component 1",x,mean,sigma) ;
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 bkg1("bkg1","Background 1",x,RooArgSet(a0,a1)) ;
46 
47  // Build exponential pdf
48  RooRealVar alpha("alpha","alpha",-1) ;
49  RooExponential bkg2("bkg2","Background 2",x,alpha) ;
50 
51  // Sum the background components into a composite background p.d.f.
52  RooRealVar bkg1frac("bkg1frac","fraction of component 1 in background",0.2,0.,1.) ;
53  RooAddPdf bkg("bkg","Signal",RooArgList(bkg1,bkg2),bkg1frac) ;
54 
55  // Sum the composite signal and background
56  RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
57  RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
58 
59 
60 
61  // Create dummy dataset that has more observables than the above pdf
62  RooRealVar y("y","y",-10,10) ;
63  RooDataSet data("data","data",RooArgSet(x,y)) ;
64 
65 
66 
67 
68  // ---------------------------------------------------
69  // B a s i c i n f o r m a t i o n r e q u e s t s
70  // ===================================================
71 
72 
73  // G e t l i s t o f o b s e r v a b l e s
74  // ---------------------------------------------
75 
76  // Get list of observables of pdf in context of a dataset
77  //
78  // Observables are define each context as the variables
79  // shared between a model and a dataset. In this case
80  // that is the variable 'x'
81 
82  RooArgSet* model_obs = model.getObservables(data) ;
83  model_obs->Print("v") ;
84 
85 
86  // G e t l i s t o f p a r a m e t e r s
87  // -------------------------------------------
88 
89  // Get list of parameters, given list of observables
90  RooArgSet* model_params = model.getParameters(x) ;
91  model_params->Print("v") ;
92 
93  // Get list of parameters, given a dataset
94  // (Gives identical results to operation above)
95  RooArgSet* model_params2 = model.getParameters(data) ;
96  model_params2->Print() ;
97 
98 
99  // G e t l i s t o f c o m p o n e n t s
100  // -------------------------------------------
101 
102  // Get list of component objects, including top-level node
103  RooArgSet* model_comps = model.getComponents() ;
104  model_comps->Print("v") ;
105 
106 
107  // -------------------------------------------------------------------------------
108  // M o d i f i c a t i o n s t o s t r u c t u r e o f c o m p o s i t e s
109  // ===============================================================================
110 
111 
112  // Create a second Gaussian
113  RooRealVar sigma2("sigma2","width of gaussians",1) ;
114  RooGaussian sig2("sig2","Signal component 1",x,mean,sigma2) ;
115 
116  // Create a sum of the original Gaussian plus the new second Gaussian
117  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
118  RooAddPdf sigsum("sigsum","sig+sig2",RooArgList(sig,sig2),sig1frac) ;
119 
120  // Construct a customizer utility to customize model
121  RooCustomizer cust(model,"cust") ;
122 
123  // Instruct the customizer to replace node 'sig' with node 'sigsum'
124  cust.replaceArg(sig,sigsum) ;
125 
126  // Build a clone of the input pdf according to the above customization
127  // instructions. Each node that requires modified is clone so that the
128  // original pdf remained untouched. The name of each cloned node is that
129  // of the original node suffixed by the name of the customizer object
130  //
131  // The returned head node own all nodes that were cloned as part of
132  // the build process so when cust_clone is deleted so will all other
133  // nodes that were created in the process.
134  RooAbsPdf* cust_clone = (RooAbsPdf*) cust.build(kTRUE) ;
135 
136  // Print structure of clone of model with sig->sigsum replacement.
137  cust_clone->Print("t") ;
138 
139 
140  delete cust_clone ;
141 
142 }
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
Exponential p.d.f.
const Double_t sigma
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
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
Double_t y[n]
Definition: legend1.C:17
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const Bool_t kTRUE
Definition: RtypesCore.h:91