Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf401_importttreethx.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit_main
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_code
9/// \macro_output
10///
11/// \date July 2008
12/// \author Wouter Verkerke
13
14#include "RooRealVar.h"
15#include "RooDataSet.h"
16#include "RooDataHist.h"
17#include "RooCategory.h"
18#include "RooGaussian.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22#include "TH1.h"
23#include "TTree.h"
24#include "TRandom.h"
25#include <map>
26
27using namespace RooFit;
28
29TH1 *makeTH1(TRandom &trnd, const char *name, double mean, double sigma);
31
32void rf401_importttreethx()
33{
34 TRandom3 trnd{};
35
36 // 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
37 // --------------------------------------------------------------------------
38
39 // Create thee ROOT TH1 histograms
40 TH1 *hh_1 = makeTH1(trnd, "hh1", 0, 3);
41 TH1 *hh_2 = makeTH1(trnd, "hh2", -3, 1);
42 TH1 *hh_3 = makeTH1(trnd, "hh3", +3, 4);
43
44 // Declare observable x
45 RooRealVar x("x", "x", -10, 10);
46
47 // Create category observable c that serves as index for the ROOT histograms
48 RooCategory c("c", "c", {{"SampleA",0}, {"SampleB",1}, {"SampleC",2}});
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 std::map<std::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
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 std::unique_ptr<RooAbsData> dsA{ds2.reduce({x, y}, "z<-5")};
103 std::unique_ptr<RooAbsData> dsB{ds2.reduce({x, y}, "abs(z)<5")};
104 std::unique_ptr<RooAbsData> dsC{ds2.reduce({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{"dsABC", "dsABC", RooArgSet(x, y), Index(c), Import("SampleA", *dsA),
108 Import("SampleB", *dsB), Import("SampleC", *dsC)};
109
110 dsABC.Print();
111}
112
113TH1 *makeTH1(TRandom &trnd, const char *name, double mean, double 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(trnd.Gaus(mean, sigma));
120 }
121 return hh;
122}
123
125{
126 // Create ROOT TTree filled with a Gaussian distribution in x and a uniform distribution in y
127
129 double *px = new double;
130 double *py = new double;
131 double *pz = new double;
132 Int_t *pi = new Int_t;
137 for (int i = 0; i < 100; i++) {
138 *px = trnd.Gaus(0, 3);
139 *py = trnd.Uniform() * 30 - 15;
140 *pz = trnd.Gaus(0, 5);
141 *pi = i % 3;
143 }
145}
#define c(i)
Definition RSha256.hxx:101
int Int_t
Definition RtypesCore.h:45
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char name[80]
Definition TGX11.cxx:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Object to represent discrete states.
Definition RooCategory.h:28
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Container class to hold unbinned data.
Definition RooDataSet.h:34
Variable that can be changed from the outside.
Definition RooRealVar.h:37
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:693
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
Random number generator class based on M.
Definition TRandom3.h:27
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
A TTree represents a columnar dataset.
Definition TTree.h:79
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg Cut(const char *cutSpec)
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64