Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf405_realtocatfuncs.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Data and categories: demonstration of real-->discrete mapping functions
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 "RooCategory.h"
16 #include "RooThresholdCategory.h"
17 #include "RooBinningCategory.h"
18 #include "Roo1DTable.h"
19 #include "RooArgusBG.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 using namespace RooFit;
24 
25 void rf405_realtocatfuncs()
26 {
27 
28  // D e f i n e p d f i n x , s a m p l e d a t a s e t i n x
29  // ------------------------------------------------------------------------
30 
31  // Define a dummy PDF in x
32  RooRealVar x("x", "x", 0, 10);
33  RooArgusBG a("a", "argus(x)", x, RooConst(10), RooConst(-1));
34 
35  // Generate a dummy dataset
36  RooDataSet *data = a.generate(x, 10000);
37 
38  // C r e a t e a t h r e s h o l d r e a l - > c a t f u n c t i o n
39  // --------------------------------------------------------------------------
40 
41  // A RooThresholdCategory is a category function that maps regions in a real-valued
42  // input observable observables to state names. At construction time a 'default'
43  // state name must be specified to which all values of x are mapped that are not
44  // otherwise assigned
45  RooThresholdCategory xRegion("xRegion", "region of x", x, "Background");
46 
47  // Specify thresholds and state assignments one-by-one.
48  // Each statement specifies that all values _below_ the given value
49  // (and above any lower specified threshold) are mapped to the
50  // category state with the given name
51  //
52  // Background | SideBand | Signal | SideBand | Background
53  // 4.23 5.23 8.23 9.23
54  xRegion.addThreshold(4.23, "Background");
55  xRegion.addThreshold(5.23, "SideBand");
56  xRegion.addThreshold(8.23, "Signal");
57  xRegion.addThreshold(9.23, "SideBand");
58 
59  // U s e t h r e s h o l d f u n c t i o n t o p l o t d a t a r e g i o n s
60  // -------------------------------------------------------------------------------------
61 
62  // Add values of threshold function to dataset so that it can be used as observable
63  data->addColumn(xRegion);
64 
65  // Make plot of data in x
66  RooPlot *xframe = x.frame(Title("Demo of threshold and binning mapping functions"));
67  data->plotOn(xframe);
68 
69  // Use calculated category to select sideband data
70  data->plotOn(xframe, Cut("xRegion==xRegion::SideBand"), MarkerColor(kRed), LineColor(kRed));
71 
72  // C r e a t e a b i n n i n g r e a l - > c a t f u n c t i o n
73  // ----------------------------------------------------------------------
74 
75  // A RooBinningCategory is a category function that maps bins of a (named) binning definition
76  // in a real-valued input observable observables to state names. The state names are automatically
77  // constructed from the variable name, the binning name and the bin number. If no binning name
78  // is specified the default binning is mapped
79 
80  x.setBins(10, "coarse");
81  RooBinningCategory xBins("xBins", "coarse bins in x", x, "coarse");
82 
83  // U s e b i n n i n g f u n c t i o n f o r t a b u l a t i o n a n d p l o t t i n g
84  // -----------------------------------------------------------------------------------------------
85 
86  // Print table of xBins state multiplicity. Note that xBins does not need to be an observable in data
87  // it can be a function of observables in data as well
88  Roo1DTable *xbtable = data->table(xBins);
89  xbtable->Print("v");
90 
91  // Add values of xBins function to dataset so that it can be used as observable
92  RooCategory *xb = (RooCategory *)data->addColumn(xBins);
93 
94  // Define range "alt" as including bins 1,3,5,7,9
95  xb->setRange("alt", "x_coarse_bin1,x_coarse_bin3,x_coarse_bin5,x_coarse_bin7,x_coarse_bin9");
96 
97  // Construct subset of data matching range "alt" but only for the first 5000 events and plot it on the frame
98  RooDataSet *dataSel = (RooDataSet *)data->reduce(CutRange("alt"), EventRange(0, 5000));
99  dataSel->plotOn(xframe, MarkerColor(kGreen), LineColor(kGreen));
100 
101  new TCanvas("rf405_realtocatfuncs", "rf405_realtocatfuncs", 600, 600);
102  xframe->SetMinimum(0.01);
103  gPad->SetLeftMargin(0.15);
104  xframe->GetYaxis()->SetTitleOffset(1.4);
105  xframe->Draw();
106 }