Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf401_importttreethx.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// Data and categories: advanced options for importing data from ROOT TTree and THx histograms
5 ///
6 /// Basic import options are demonstrated in rf102_dataimport.C
7 ///
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 #include "RooRealVar.h"
13 #include "RooDataSet.h"
14 #include "RooDataHist.h"
15 #include "RooCategory.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 #include "TH1.h"
22 #include "TTree.h"
23 #include "TRandom.h"
24 #include <map>
25 
26 using namespace RooFit;
27 
28 TH1 *makeTH1(const char *name, Double_t mean, Double_t sigma);
29 TTree *makeTTree();
30 
31 void rf401_importttreethx()
32 {
33  // I m p o r t m u l t i p l e T H 1 i n t o a R o o D a t a H i s t
34  // --------------------------------------------------------------------------
35 
36  // Create thee ROOT TH1 histograms
37  TH1 *hh_1 = makeTH1("hh1", 0, 3);
38  TH1 *hh_2 = makeTH1("hh2", -3, 1);
39  TH1 *hh_3 = makeTH1("hh3", +3, 4);
40 
41  // Declare observable x
42  RooRealVar x("x", "x", -10, 10);
43 
44  // Create category observable c that serves as index for the ROOT histograms
45  RooCategory c("c", "c");
46  c.defineType("SampleA");
47  c.defineType("SampleB");
48  c.defineType("SampleC");
49 
50  // Create a binned dataset that imports contents of all TH1 mapped by index category c
51  RooDataHist *dh = new RooDataHist("dh", "dh", x, Index(c), Import("SampleA", *hh_1), Import("SampleB", *hh_2),
52  Import("SampleC", *hh_3));
53  dh->Print();
54 
55  // Alternative constructor form for importing multiple histograms
56  map<string, TH1 *> hmap;
57  hmap["SampleA"] = hh_1;
58  hmap["SampleB"] = hh_2;
59  hmap["SampleC"] = hh_3;
60  RooDataHist *dh2 = new RooDataHist("dh", "dh", x, c, hmap);
61  dh2->Print();
62 
63  // I m p o r t i n g a T T r e e i n t o a R o o D a t a S e t w i t h c u t s
64  // -----------------------------------------------------------------------------------------
65 
66  TTree *tree = makeTTree();
67 
68  // Define observables y,z
69  RooRealVar y("y", "y", -10, 10);
70  RooRealVar z("z", "z", -10, 10);
71 
72  // Import only observables (y,z)
73  RooDataSet ds("ds", "ds", RooArgSet(x, y), Import(*tree));
74  ds.Print();
75 
76  // Import observables (x,y,z) but only event for which (y+z<0) is true
77  RooDataSet ds2("ds2", "ds2", RooArgSet(x, y, z), Import(*tree), Cut("y+z<0"));
78  ds2.Print();
79 
80  // I m p o r t i n g i n t e g e r T T r e e b r a n c h e s
81  // ---------------------------------------------------------------
82 
83  // Import integer tree branch as RooRealVar
84  RooRealVar i("i", "i", 0, 5);
85  RooDataSet ds3("ds3", "ds3", RooArgSet(i, x), Import(*tree));
86  ds3.Print();
87 
88  // Define category i
89  RooCategory icat("i", "i");
90  icat.defineType("State0", 0);
91  icat.defineType("State1", 1);
92 
93  // Import integer tree branch as RooCategory (only events with i==0 and i==1
94  // will be imported as those are the only defined states)
95  RooDataSet ds4("ds4", "ds4", RooArgSet(icat, x), Import(*tree));
96  ds4.Print();
97 
98  // I m p o r t m u l t i p l e R o o D a t a S e t s i n t o a R o o D a t a S e t
99  // ----------------------------------------------------------------------------------------
100 
101  // Create three RooDataSets in (y,z)
102  RooDataSet *dsA = (RooDataSet *)ds2.reduce(RooArgSet(x, y), "z<-5");
103  RooDataSet *dsB = (RooDataSet *)ds2.reduce(RooArgSet(x, y), "abs(z)<5");
104  RooDataSet *dsC = (RooDataSet *)ds2.reduce(RooArgSet(x, y), "z>5");
105 
106  // Create a dataset that imports contents of all the above datasets mapped by index category c
107  RooDataSet *dsABC = new RooDataSet("dsABC", "dsABC", RooArgSet(x, y), Index(c), Import("SampleA", *dsA),
108  Import("SampleB", *dsB), Import("SampleC", *dsC));
109 
110  dsABC->Print();
111 }
112 
113 TH1 *makeTH1(const char *name, Double_t mean, Double_t sigma)
114 {
115  // Create ROOT TH1 filled with a Gaussian distribution
116 
117  TH1D *hh = new TH1D(name, name, 100, -10, 10);
118  for (int i = 0; i < 1000; i++) {
119  hh->Fill(gRandom->Gaus(mean, sigma));
120  }
121  return hh;
122 }
123 
124 TTree *makeTTree()
125 {
126  // Create ROOT TTree filled with a Gaussian distribution in x and a uniform distribution in y
127 
128  TTree *tree = new TTree("tree", "tree");
129  Double_t *px = new Double_t;
130  Double_t *py = new Double_t;
131  Double_t *pz = new Double_t;
132  Int_t *pi = new Int_t;
133  tree->Branch("x", px, "x/D");
134  tree->Branch("y", py, "y/D");
135  tree->Branch("z", pz, "z/D");
136  tree->Branch("i", pi, "i/I");
137  for (int i = 0; i < 100; i++) {
138  *px = gRandom->Gaus(0, 3);
139  *py = gRandom->Uniform() * 30 - 15;
140  *pz = gRandom->Gaus(0, 5);
141  *pi = i % 3;
142  tree->Fill();
143  }
144  return tree;
145 }