Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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_code
7/// \macro_output
8///
9/// \date July 2008
10/// \author Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooDataSet.h"
14#include "RooGaussian.h"
15#include "RooProdPdf.h"
16#include "RooAddPdf.h"
17#include "RooPolynomial.h"
18#include "TCanvas.h"
19#include "TAxis.h"
20#include "RooPlot.h"
21#include "RooFitResult.h"
22using namespace RooFit;
23
25{
26
27 // C r e a t e 2 D p d f a n d d a t a
28 // -------------------------------------------
29
30 // Define observables x,y
31 RooRealVar x("x", "x", -10, 10);
32 RooRealVar y("y", "y", -10, 10);
33
34 // Construct the signal pdf gauss(x)*gauss(y)
35 RooRealVar mx("mx", "mx", 1, -10, 10);
36 RooRealVar my("my", "my", 1, -10, 10);
37
38 RooGaussian gx("gx", "gx", x, mx, 1.0);
39 RooGaussian gy("gy", "gy", y, my, 1.0);
40
41 RooProdPdf sig("sig", "sig", gx, gy);
42
43 // Construct the background pdf (flat in x,y)
44 RooPolynomial px("px", "px", x);
45 RooPolynomial py("py", "py", y);
46 RooProdPdf bkg("bkg", "bkg", px, py);
47
48 // Construct the composite model sig+bkg
49 RooRealVar f("f", "f", 0., 1.);
50 RooAddPdf model("model", "model", RooArgList(sig, bkg), f);
51
52 // Sample 10000 events in (x,y) from the model
53 std::unique_ptr<RooDataSet> modelData{model.generate({x, y}, 10000)};
54
55 // 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
56 // -------------------------------------------------------------------
57
58 // Construct the SideBand1,SideBand2,Signal regions
59 //
60 // |
61 // +-------------+-----------+
62 // | | |
63 // | Side | Sig |
64 // | Band1 | nal |
65 // | | |
66 // --+-------------+-----------+--
67 // | |
68 // | Side |
69 // | Band2 |
70 // | |
71 // +-------------+-----------+
72 // |
73
74 x.setRange("SB1", -10, +10);
75 y.setRange("SB1", -10, 0);
76
77 x.setRange("SB2", -10, 0);
78 y.setRange("SB2", 0, +10);
79
80 x.setRange("SIG", 0, +10);
81 y.setRange("SIG", 0, +10);
82
83 x.setRange("FULL", -10, +10);
84 y.setRange("FULL", -10, +10);
85
86 // 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
87 // -------------------------------------------------------------------------------------
88
89 // Perform fit in SideBand1 region (RooAddPdf coefficients will be interpreted in full range)
90 std::unique_ptr<RooFitResult> r_sb1{model.fitTo(*modelData, Range("SB1"), Save(), PrintLevel(-1))};
91
92 // Perform fit in SideBand2 region (RooAddPdf coefficients will be interpreted in full range)
93 std::unique_ptr<RooFitResult> r_sb2{model.fitTo(*modelData, Range("SB2"), Save(), PrintLevel(-1))};
94
95 // 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
96 // -----------------------------------------------------------------------------
97
98 // Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
99 // (RooAddPdf coefficients will be interpreted in full range)
100 std::unique_ptr<RooFitResult> r_sb12{model.fitTo(*modelData, Range("SB1,SB2"), Save(), PrintLevel(-1))};
101
102 // Print results for comparison
103 r_sb1->Print();
104 r_sb2->Print();
105 r_sb12->Print();
106}
#define f(i)
Definition RSha256.hxx:104
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char mx
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooPolynomial implements a polynomial p.d.f of the form.
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
Ta Range(0, 0, 1, 1)