Logo ROOT   6.10/09
Reference Guide
rf312_multirangefit.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #312
4 //
5 // Performing fits in multiple (disjoint) ranges in one or more dimensions
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooProdPdf.h"
20 #include "RooAddPdf.h"
21 #include "RooPolynomial.h"
22 #include "TCanvas.h"
23 #include "RooPlot.h"
24 #include "RooFitResult.h"
25 using namespace RooFit ;
26 
27 
28 class TestBasic312 : public RooFitTestUnit
29 {
30 public:
31  TestBasic312(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fit in multiple rectangular ranges",refFile,writeRef,verbose) {} ;
32  Bool_t testCode() {
33 
34  // C r e a t e 2 D p d f a n d d a t a
35  // -------------------------------------------
36 
37  // Define observables x,y
38  RooRealVar x("x","x",-10,10) ;
39  RooRealVar y("y","y",-10,10) ;
40 
41  // Construct the signal pdf gauss(x)*gauss(y)
42  RooRealVar mx("mx","mx",1,-10,10) ;
43  RooRealVar my("my","my",1,-10,10) ;
44 
45  RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;
46  RooGaussian gy("gy","gy",y,my,RooConst(1)) ;
47 
48  RooProdPdf sig("sig","sig",gx,gy) ;
49 
50  // Construct the background pdf (flat in x,y)
51  RooPolynomial px("px","px",x) ;
52  RooPolynomial py("py","py",y) ;
53  RooProdPdf bkg("bkg","bkg",px,py) ;
54 
55  // Construct the composite model sig+bkg
56  RooRealVar f("f","f",0.,1.) ;
57  RooAddPdf model("model","model",RooArgList(sig,bkg),f) ;
58 
59  // Sample 10000 events in (x,y) from the model
60  RooDataSet* modelData = model.generate(RooArgSet(x,y),10000) ;
61 
62 
63 
64  // 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
65  // -------------------------------------------------------------------
66 
67  // Construct the SideBand1,SideBand2,Signal regions
68  //
69  // |
70  // +-------------+-----------+
71  // | | |
72  // | Side | Sig |
73  // | Band1 | nal |
74  // | | |
75  // --+-------------+-----------+--
76  // | |
77  // | Side |
78  // | Band2 |
79  // | |
80  // +-------------+-----------+
81  // |
82 
83  x.setRange("SB1",-10,+10) ;
84  y.setRange("SB1",-10,0) ;
85 
86  x.setRange("SB2",-10,0) ;
87  y.setRange("SB2",0,+10) ;
88 
89  x.setRange("SIG",0,+10) ;
90  y.setRange("SIG",0,+10) ;
91 
92  x.setRange("FULL",-10,+10) ;
93  y.setRange("FULL",-10,+10) ;
94 
95 
96  // 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
97  // -------------------------------------------------------------------------------------
98 
99  // Perform fit in SideBand1 region (RooAddPdf coefficients will be interpreted in full range)
100  RooFitResult* r_sb1 = model.fitTo(*modelData,Range("SB1"),Save()) ;
101 
102  // Perform fit in SideBand2 region (RooAddPdf coefficients will be interpreted in full range)
103  RooFitResult* r_sb2 = model.fitTo(*modelData,Range("SB2"),Save()) ;
104 
105 
106 
107  // 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
108  // -----------------------------------------------------------------------------
109 
110  // Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
111  // (RooAddPdf coefficients will be interpreted in full range)
112  RooFitResult* r_sb12 = model.fitTo(*modelData,Range("SB1,SB2"),Save()) ;
113 
114 
115  regResult(r_sb1,"rf312_fit_sb1") ;
116  regResult(r_sb2,"rf312_fit_sb2") ;
117  regResult(r_sb12,"rf312_fit_sb12") ;
118 
119  delete modelData ;
120 
121  return kTRUE ;
122 
123  }
124 } ;
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
double f(double x)
Double_t y[n]
Definition: legend1.C:17
RooCmdArg Save(Bool_t flag=kTRUE)
RooConstVar & RooConst(Double_t val)
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
const Bool_t kTRUE
Definition: RtypesCore.h:91