Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
multidimSampling.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// \notebook
4/// Example of random number generation by sampling a multi-dim distribution using the DistSampler class.
5///
6/// NOTE: This tutorial must be run with ACLIC
7///
8/// \macro_image
9/// \macro_code
10///
11/// \author Lorenzo Moneta
12
13// function (a 4d gaussian)
14#include "TMath.h"
15#include "TF2.h"
16#include "TStopwatch.h"
17#include "Math/DistSampler.h"
20#include "Math/Factory.h"
21
22#include "TKDTreeBinning.h"
23
24#include "TTree.h"
25#include "TFile.h"
26#include "TMatrixDSym.h"
27#include "TVectorD.h"
28#include "TCanvas.h"
29#include <cmath>
30
31// Gauss ND function
32// make a class in order to avoid constructing the
33// matrices for every call
34// This however requires that the code must be compiled with ACLIC
35
36bool debug = false;
37
38// Define the GausND structure
39struct GausND {
40
41 TVectorD X;
44
45 GausND( int dim ) :
46 X(TVectorD(dim)),
47 Mu(TVectorD(dim)),
48 CovMat(TMatrixDSym(dim) )
49 {}
50
51 double operator() (double *x, double *p) {
52 // 4 parameters
53 int dim = X.GetNrows();
54 int k = 0;
55 for (int i = 0; i<dim; ++i) { X[i] = x[i] - p[k]; k++; }
56 for (int i = 0; i<dim; ++i) {
57 CovMat(i,i) = p[k]*p[k];
58 k++;
59 }
60 for (int i = 0; i<dim; ++i) {
61 for (int j = i+1; j<dim; ++j) {
62 // p now are the correlations N(N-1)/2
63 CovMat(i,j) = p[k]*sqrt(CovMat(i,i)*CovMat(j,j));
64 CovMat(j,i) = CovMat(i,j);
65 k++;
66 }
67 }
68 if (debug) {
69 X.Print();
70 CovMat.Print();
71 }
72
73 double det = CovMat.Determinant();
74 if (det <= 0) {
75 Fatal("GausND","Determinant is <= 0 det = %f",det);
76 CovMat.Print();
77 return 0;
78 }
79 double norm = std::pow( 2. * TMath::Pi(), dim/2) * sqrt(det);
80 // compute the gaussians
81 CovMat.Invert();
82 double fval = std::exp( - 0.5 * CovMat.Similarity(X) )/ norm;
83
84 if (debug) {
85 std::cout << "det " << det << std::endl;
86 std::cout << "norm " << norm << std::endl;
87 std::cout << "fval " << fval << std::endl;
88 }
89
90 return fval;
91 }
92};
93
94// Use the Math namespace
95using namespace ROOT::Math;
96
97void multidimSampling() {
98
99
100 const int N = 10000;
101 /*const int NBin = 1000;*/
102 const int DIM = 4;
103
104 double xmin[] = {-10,-10,-10, -10};
105 double xmax[] = { 10, 10, 10, 10};
106 double par0[] = { 1., -1., 2, 0, // the gaussian mu
107 1, 2, 1, 3, // the sigma
108 0.5,0.,0.,0.,0.,0.8 }; // the correlation
109
110 const int NPAR = DIM + DIM*(DIM+1)/2; // 14 in the 4d case
111 // generate the sample
112 GausND gaus4d(4);
113 TF1 * f = new TF1("functionND",gaus4d,0,1,14);
114 f->SetParameters(par0);
115
116 double x0[] = {0,0,0,0};
117 // for debugging
118 if (debug) f->EvalPar(x0,nullptr);
119 debug = false;
120
122 for (int i = 0; i < NPAR; ++i ) {
123 if (i < DIM) f->SetParName(i, name.Format("mu_%d",i+1) );
124 else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-DIM+1) );
125 else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-2*DIM+1) );
126 }
127
128 /*ROOT::Math::DistSamplerOptions::SetDefaultSampler("Foam");*/
129 DistSampler * sampler = Factory::CreateDistSampler();
130 if (sampler == nullptr) {
131 Info("multidimSampling","Default sampler %s is not available try with Foam ",
134 }
135 sampler = Factory::CreateDistSampler();
136 if (sampler == nullptr) {
137 Error("multidimSampling","Foam sampler is not available - exit ");
138 return;
139 }
140
141 sampler->SetFunction(*f,DIM);
142 sampler->SetRange(xmin,xmax);
143 bool ret = sampler->Init();
144
145 std::vector<double> data1(DIM*N);
146 double v[DIM];
148
149 if (!ret) {
150 Error("Sampler::Init","Error initializing unuran sampler");
151 return;
152 }
153
154 // generate the data
155 w.Start();
156 for (int i = 0; i < N; ++i) {
157 sampler->Sample(v);
158 for (int j = 0; j < DIM; ++j)
159 data1[N*j + i] = v[j];
160 }
161 w.Stop();
162 w.Print();
163
164 // fill tree with data
165 TFile * file = new TFile("multiDimSampling.root","RECREATE");
166 double x[DIM];
167 TTree * t1 = new TTree("t1","Tree from Unuran");
168 t1->Branch("x",x,"x[4]/D");
169 for (int i = 0; i < N; ++i) {
170 for (int j = 0; j < DIM; ++j) {
171 x[j] = data1[i+N*j];
172 }
173 t1->Fill();
174 }
175
176 // plot the data
177 t1->Draw("x[0]:x[1]:x[2]:x[3]","","candle");
178 TCanvas * c2 = new TCanvas();
179 c2->Divide(3,2);
180 int ic=1;
181 c2->cd(ic++);
182 t1->Draw("x[0]:x[1]");
183 c2->cd(ic++);
184 t1->Draw("x[0]:x[2]");
185 c2->cd(ic++);
186 t1->Draw("x[0]:x[3]");
187 c2->cd(ic++);
188 t1->Draw("x[1]:x[2]");
189 c2->cd(ic++);
190 t1->Draw("x[1]:x[3]");
191 c2->cd(ic++);
192 t1->Draw("x[2]:x[3]");
193
194 t1->Write();
195 file->Close();
196
197}
#define f(i)
Definition RSha256.hxx:104
#define X(type, name)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Fatal(const char *location, const char *msgfmt,...)
Use this function in case of a fatal error. It will abort the program.
Definition TError.cxx:244
#define N
winID h TVirtualViewer3D TVirtualGLPainter p
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
TRObject operator()(const T1 &t1) const
static void SetDefaultSampler(const char *type)
static const std::string & DefaultSampler()
Interface class for generic sampling of a distribution, i.e.
Definition DistSampler.h:58
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:234
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
void Close(Option_t *option="") override
Close a file.
Definition TFile.cxx:955
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:139
A TTree represents a columnar dataset.
Definition TTree.h:84
std::ostream & Info()
Definition hadd.cxx:171
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
constexpr Double_t Pi()
Definition TMath.h:37
auto * t1
Definition textangle.C:20