Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf203_ranges.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Addition and convolution: fitting and plotting in sub ranges
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooPolynomial.h"
16 #include "RooAddPdf.h"
17 #include "RooFitResult.h"
18 #include "RooPlot.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "TH1.h"
22 using namespace RooFit;
23 
24 void rf203_ranges()
25 {
26  // S e t u p m o d e l
27  // ---------------------
28 
29  // Construct observables x
30  RooRealVar x("x", "x", -10, 10);
31 
32  // Construct gaussx(x,mx,1)
33  RooRealVar mx("mx", "mx", 0, -10, 10);
34  RooGaussian gx("gx", "gx", x, mx, RooConst(1));
35 
36  // Construct px = 1 (flat in x)
37  RooPolynomial px("px", "px", x);
38 
39  // Construct model = f*gx + (1-f)px
40  RooRealVar f("f", "f", 0., 1.);
41  RooAddPdf model("model", "model", RooArgList(gx, px), f);
42 
43  // Generated 10000 events in (x,y) from p.d.f. model
44  RooDataSet *modelData = model.generate(x, 10000);
45 
46  // F i t f u l l r a n g e
47  // ---------------------------
48 
49  // Fit p.d.f to all data
50  RooFitResult *r_full = model.fitTo(*modelData, Save(kTRUE));
51 
52  // F i t p a r t i a l r a n g e
53  // ----------------------------------
54 
55  // Define "signal" range in x as [-3,3]
56  x.setRange("signal", -3, 3);
57 
58  // Fit p.d.f only to data in "signal" range
59  RooFitResult *r_sig = model.fitTo(*modelData, Save(kTRUE), Range("signal"));
60 
61  // P l o t / p r i n t r e s u l t s
62  // ---------------------------------------
63 
64  // Make plot frame in x and add data and fitted model
65  RooPlot *frame = x.frame(Title("Fitting a sub range"));
66  modelData->plotOn(frame);
67  model.plotOn(frame, Range("Full"), LineStyle(kDashed), LineColor(kRed)); // Add shape in full ranged dashed
68  model.plotOn(frame); // By default only fitted range is shown
69 
70  // Print fit results
71  cout << "result of fit on all data " << endl;
72  r_full->Print();
73  cout << "result of fit in in signal region (note increased error on signal fraction)" << endl;
74  r_sig->Print();
75 
76  // Draw frame on canvas
77  new TCanvas("rf203_ranges", "rf203_ranges", 600, 600);
78  gPad->SetLeftMargin(0.15);
79  frame->GetYaxis()->SetTitleOffset(1.4);
80  frame->Draw();
81 
82  return;
83 }