Logo ROOT   6.10/09
Reference Guide
rf203_ranges.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #203
4 //
5 // Fitting and plotting in sub ranges
6 //
7 //
8 // 07/2008 - Wouter Verkerke
9 //
10 /////////////////////////////////////////////////////////////////////////
11 
12 #ifndef __CINT__
13 #include "RooGlobalFunc.h"
14 #endif
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "RooPolynomial.h"
19 #include "RooAddPdf.h"
20 #include "RooFitResult.h"
21 #include "RooPlot.h"
22 #include "TCanvas.h"
23 #include "TH1.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic203 : public RooFitTestUnit
28 {
29 public:
30  TestBasic203(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Basic fitting and plotting in ranges",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // S e t u p m o d e l
34  // ---------------------
35 
36  // Construct observables x
37  RooRealVar x("x","x",-10,10) ;
38 
39  // Construct gaussx(x,mx,1)
40  RooRealVar mx("mx","mx",0,-10,10) ;
41  RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;
42 
43  // Construct px = 1 (flat in x)
44  RooPolynomial px("px","px",x) ;
45 
46  // Construct model = f*gx + (1-f)px
47  RooRealVar f("f","f",0.,1.) ;
48  RooAddPdf model("model","model",RooArgList(gx,px),f) ;
49 
50  // Generated 10000 events in (x,y) from p.d.f. model
51  RooDataSet* modelData = model.generate(x,10000) ;
52 
53  // F i t f u l l r a n g e
54  // ---------------------------
55 
56  // Fit p.d.f to all data
57  RooFitResult* r_full = model.fitTo(*modelData,Save(kTRUE)) ;
58 
59 
60  // F i t p a r t i a l r a n g e
61  // ----------------------------------
62 
63  // Define "signal" range in x as [-3,3]
64  x.setRange("signal",-3,3) ;
65 
66  // Fit p.d.f only to data in "signal" range
67  RooFitResult* r_sig = model.fitTo(*modelData,Save(kTRUE),Range("signal")) ;
68 
69 
70  // P l o t / p r i n t r e s u l t s
71  // ---------------------------------------
72 
73  // Make plot frame in x and add data and fitted model
74  RooPlot* frame = x.frame(Title("Fitting a sub range")) ;
75  modelData->plotOn(frame) ;
76  model.plotOn(frame,Range("Full"),LineStyle(kDashed),LineColor(kRed)) ; // Add shape in full ranged dashed
77  model.plotOn(frame) ; // By default only fitted range is shown
78 
79  regPlot(frame,"rf203_plot") ;
80  regResult(r_full,"rf203_r_full") ;
81  regResult(r_sig,"rf203_r_sig") ;
82 
83  delete modelData ;
84  return kTRUE;
85  }
86 } ;
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
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:552
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
Double_t x[n]
Definition: legend1.C:17
RooCmdArg LineStyle(Style_t style)
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
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
double f(double x)
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