ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
multidimSampling.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// Example of sampling a multi-dim distribution using the DistSampler class
4 /// NOTE: This tutorial must be run with ACLIC
5 ///
6 /// \macro_image
7 /// \macro_code
8 ///
9 /// \author Lorenzo Moneta
10 
11 // function (a 4d gaussian)
12 #include "TMath.h"
13 #include "TF2.h"
14 #include "TStopwatch.h"
15 #include "Math/DistSampler.h"
17 #include "Math/MinimizerOptions.h"
18 #include "Math/Factory.h"
19 
20 #include "TKDTreeBinning.h"
21 
22 #include "TTree.h"
23 #include "TFile.h"
24 #include "TMatrixDSym.h"
25 #include "TVectorD.h"
26 #include "TCanvas.h"
27 #include <cmath>
28 
29 // Gauss ND function
30 // make a class in order to avoid constructing the
31 // matrices for every call
32 // This however requires that the code must be compiled with ACLIC
33 
34 bool debug = false;
35 struct GausND {
36 
37  TVectorD X;
38  TVectorD Mu;
39  TMatrixDSym CovMat;
40 
41  GausND( int dim ) :
42  X(TVectorD(dim)),
43  Mu(TVectorD(dim)),
44  CovMat(TMatrixDSym(dim) )
45  {}
46 
47  double operator() (double *x, double *p) {
48  // 4 parameters
49  int dim = X.GetNrows();
50  int k = 0;
51  for (int i = 0; i<dim; ++i) { X[i] = x[i] - p[k]; k++; }
52  for (int i = 0; i<dim; ++i) {
53  CovMat(i,i) = p[k]*p[k];
54  k++;
55  }
56  for (int i = 0; i<dim; ++i) {
57  for (int j = i+1; j<dim; ++j) {
58  // p now are the correlations N(N-1)/2
59  CovMat(i,j) = p[k]*sqrt(CovMat(i,i)*CovMat(j,j));
60  CovMat(j,i) = CovMat(i,j);
61  k++;
62  }
63  }
64  if (debug) {
65  X.Print();
66  CovMat.Print();
67  }
68 
69  double det = CovMat.Determinant();
70  if (det <= 0) {
71  Fatal("GausND","Determinant is <= 0 det = %f",det);
72  CovMat.Print();
73  return 0;
74  }
75  double norm = std::pow( 2. * TMath::Pi(), dim/2) * sqrt(det);
76  // compute the gaussians
77  CovMat.Invert();
78  double fval = std::exp( - 0.5 * CovMat.Similarity(X) )/ norm;
79 
80  if (debug) {
81  std::cout << "det " << det << std::endl;
82  std::cout << "norm " << norm << std::endl;
83  std::cout << "fval " << fval << std::endl;
84  }
85 
86  return fval;
87  }
88 };
89 
90 using namespace ROOT::Math;
91 
92 void multidimSampling() {
93 
94 #ifdef __CINT__
95  std::cout << "DO NOT RUN WITH CINT:" << std::endl;
96  std::cout << "we are using a custom function which requires" << std::endl;
97  std::cout << "that this tutorial must be compiled with ACLIC" << std::endl;
98  return;
99 #endif
100 
101  const int N = 100000;
102  //const int NBin = 1000;
103  const int DIM = 4;
104 
105  double xmin[] = {-10,-10,-10, -10};
106  double xmax[] = { 10, 10, 10, 10};
107  double par0[] = { 1., -1., 2, 0, // the gaussian mu
108  1, 2, 1, 3, // the sigma
109  0.5,0.,0.,0.,0.,0.8 }; // the correlation
110 
111  const int NPAR = DIM + DIM*(DIM+1)/2; // 14 in the 4d case
112  // generate the sample
113  GausND gaus4d(4);
114  TF1 * f = new TF1("functionND",gaus4d,0,1,14);
115  f->SetParameters(par0);
116 
117  double x0[] = {0,0,0,0};
118  // for debugging
119  if (debug) f->EvalPar(x0,0);
120  debug = false;
121 
122  TString name;
123  for (int i = 0; i < NPAR; ++i ) {
124  if (i < DIM) f->SetParName(i, name.Format("mu_%d",i+1) );
125  else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-DIM+1) );
126  else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-2*DIM+1) );
127  }
128 
129  //ROOT::Math::DistSamplerOptions::SetDefaultSampler("Foam");
131  if (sampler == 0) {
132  Info("multidimSampling","Default sampler %s is not available try with Foam ",
135  }
136  sampler = Factory::CreateDistSampler();
137  if (sampler == 0) {
138  Error("multidimSampling","Foam sampler is not available - exit ");
139  return;
140  }
141 
142  sampler->SetFunction(*f,DIM);
143  sampler->SetRange(xmin,xmax);
144  bool ret = sampler->Init();
145 
146  std::vector<double> data1(DIM*N);
147  double v[DIM];
148  TStopwatch w;
149 
150  if (!ret) {
151  Error("Sampler::Init","Error initializing unuran sampler");
152  return;
153  }
154 
155  // generate the data
156  w.Start();
157  for (int i = 0; i < N; ++i) {
158  sampler->Sample(v);
159  for (int j = 0; j < DIM; ++j)
160  data1[N*j + i] = v[j];
161  }
162  w.Stop();
163  w.Print();
164 
165  // fill tree with data
166  TFile * file = new TFile("multiDimSampling.root","RECREATE");
167  double x[DIM];
168  TTree * t1 = new TTree("t1","Tree from Unuran");
169  t1->Branch("x",x,"x[4]/D");
170  for (int i = 0; i < N; ++i) {
171  for (int j = 0; j < DIM; ++j) {
172  x[j] = data1[i+N*j];
173  }
174  t1->Fill();
175  }
176 
177  // plot the data
178  t1->Draw("x[0]:x[1]:x[2]:x[3]","","candle");
179  TCanvas * c2 = new TCanvas();
180  c2->Divide(3,2);
181  int ic=1;
182  c2->cd(ic++);
183  t1->Draw("x[0]:x[1]");
184  c2->cd(ic++);
185  t1->Draw("x[0]:x[2]");
186  c2->cd(ic++);
187  t1->Draw("x[0]:x[3]");
188  c2->cd(ic++);
189  t1->Draw("x[1]:x[2]");
190  c2->cd(ic++);
191  t1->Draw("x[1]:x[3]");
192  c2->cd(ic++);
193  t1->Draw("x[2]:x[3]");
194 
195  t1->Write();
196  file->Close();
197 
198 }
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:217
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:56
void Fatal(const char *location, const char *msgfmt,...)
static void SetDefaultSampler(const char *type)
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4306
const double * Sample()
sample one event and rerturning array x with coordinates
Definition: DistSampler.h:173
#define N
Basic string class.
Definition: TString.h:137
double sqrt(double)
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:75
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2321
void SetFunction(Function &func, unsigned int dim)
set the parent function distribution to use for sampling (generic case)
Definition: DistSampler.h:76
double pow(double, double)
void Info(const char *location, const char *msgfmt,...)
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TTree.cxx:8780
tuple w
Definition: qtexample.py:51
static ROOT::Math::DistSampler * CreateDistSampler(const std::string &samplerType="")
static method to create the distribution sampler class given a string specifying the type Supported s...
Definition: Factory.cxx:167
Double_t Pi()
Definition: TMath.h:44
void SetRange(double xmin, double xmax, int icoord=0)
set range in a given dimension
Definition: DistSampler.cxx:40
TRObject operator()(const T1 &t1) const
Interface class for generic sampling of a distribution, i.e.
Definition: DistSampler.h:61
static const std::string & DefaultSampler()
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:360
#define name(a, b)
Definition: linkTestLib0.cpp:5
virtual Int_t Branch(TCollection *list, Int_t bufsize=32000, Int_t splitlevel=99, const char *name="")
Create one branch for each element in the collection.
Definition: TTree.cxx:1623
virtual bool Init(const char *="")
initialize the generators with the given algorithm Implemented by derived classes who needs it (like ...
Definition: DistSampler.h:105
bool debug
A TTree object has a header with a name and a title.
Definition: TTree.h:98
const int NPAR
Definition: ErrorIntegral.C:22
double exp(double)
void Error(ErrorHandler_t func, int code, const char *va_(fmt),...)
Write error message and call a handler, if required.
Stopwatch class.
Definition: TStopwatch.h:30