Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf312_multirangefit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// Multidimensional models: performing fits in multiple (disjoint) ranges in one or more dimensions
5 ///
6 /// \macro_output
7 /// \macro_code
8 /// \author 07/2008 - Wouter Verkerke
9 
10 #include "RooRealVar.h"
11 #include "RooDataSet.h"
12 #include "RooGaussian.h"
13 #include "RooConstVar.h"
14 #include "RooProdPdf.h"
15 #include "RooAddPdf.h"
16 #include "RooPolynomial.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "RooPlot.h"
20 #include "RooFitResult.h"
21 using namespace RooFit;
22 
23 void rf312_multirangefit()
24 {
25 
26  // C r e a t e 2 D p d f a n d d a t a
27  // -------------------------------------------
28 
29  // Define observables x,y
30  RooRealVar x("x", "x", -10, 10);
31  RooRealVar y("y", "y", -10, 10);
32 
33  // Construct the signal pdf gauss(x)*gauss(y)
34  RooRealVar mx("mx", "mx", 1, -10, 10);
35  RooRealVar my("my", "my", 1, -10, 10);
36 
37  RooGaussian gx("gx", "gx", x, mx, RooConst(1));
38  RooGaussian gy("gy", "gy", y, my, RooConst(1));
39 
40  RooProdPdf sig("sig", "sig", gx, gy);
41 
42  // Construct the background pdf (flat in x,y)
43  RooPolynomial px("px", "px", x);
44  RooPolynomial py("py", "py", y);
45  RooProdPdf bkg("bkg", "bkg", px, py);
46 
47  // Construct the composite model sig+bkg
48  RooRealVar f("f", "f", 0., 1.);
49  RooAddPdf model("model", "model", RooArgList(sig, bkg), f);
50 
51  // Sample 10000 events in (x,y) from the model
52  RooDataSet *modelData = model.generate(RooArgSet(x, y), 10000);
53 
54  // D e f i n e s i g n a l a n d s i d e b a n d r e g i o n s
55  // -------------------------------------------------------------------
56 
57  // Construct the SideBand1,SideBand2,Signal regions
58  //
59  // |
60  // +-------------+-----------+
61  // | | |
62  // | Side | Sig |
63  // | Band1 | nal |
64  // | | |
65  // --+-------------+-----------+--
66  // | |
67  // | Side |
68  // | Band2 |
69  // | |
70  // +-------------+-----------+
71  // |
72 
73  x.setRange("SB1", -10, +10);
74  y.setRange("SB1", -10, 0);
75 
76  x.setRange("SB2", -10, 0);
77  y.setRange("SB2", 0, +10);
78 
79  x.setRange("SIG", 0, +10);
80  y.setRange("SIG", 0, +10);
81 
82  x.setRange("FULL", -10, +10);
83  y.setRange("FULL", -10, +10);
84 
85  // P e r f o r m f i t s i n i n d i v i d u a l s i d e b a n d r e g i o n s
86  // -------------------------------------------------------------------------------------
87 
88  // Perform fit in SideBand1 region (RooAddPdf coefficients will be interpreted in full range)
89  RooFitResult *r_sb1 = model.fitTo(*modelData, Range("SB1"), Save());
90 
91  // Perform fit in SideBand2 region (RooAddPdf coefficients will be interpreted in full range)
92  RooFitResult *r_sb2 = model.fitTo(*modelData, Range("SB2"), Save());
93 
94  // P e r f o r m f i t s i n j o i n t s i d e b a n d r e g i o n s
95  // -----------------------------------------------------------------------------
96 
97  // Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
98  // (RooAddPdf coefficients will be interpreted in full range)
99  RooFitResult *r_sb12 = model.fitTo(*modelData, Range("SB1,SB2"), Save());
100 
101  // Print results for comparison
102  r_sb1->Print();
103  r_sb2->Print();
104  r_sb12->Print();
105 }