00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "TStopwatch.h"
00010 #include "TUnuran.h"
00011 #include "TUnuranMultiContDist.h"
00012
00013 #include "TH3.h"
00014 #include "TF3.h"
00015 #include "TCanvas.h"
00016 #include "TMath.h"
00017
00018 #include "TRandom3.h"
00019 #include "TSystem.h"
00020 #include "TApplication.h"
00021 #include "TError.h"
00022
00023
00024
00025
00026
00027 #include <iostream>
00028 #include <cmath>
00029 #include <cassert>
00030
00031 #include <vector>
00032
00033 using std::cout;
00034 using std::endl;
00035
00036 int n;
00037
00038 bool useRandomSeed = false;
00039
00040 double gaus3d(double *x, double *p) {
00041
00042 double sigma_x = p[0];
00043 double sigma_y = p[1];
00044 double sigma_z = p[2];
00045 double rho = p[3];
00046 double u = x[0] / sigma_x ;
00047 double v = x[1] / sigma_y ;
00048 double w = x[2] / sigma_z ;
00049 double c = 1 - rho*rho ;
00050 double result = (1 / (2 * TMath::Pi() * sigma_x * sigma_y * sigma_z * sqrt(c)))
00051 * exp (-(u * u - 2 * rho * u * v + v * v + w*w) / (2 * c));
00052 return result;
00053 }
00054 double log_gaus3d(double *x, double *p) {
00055 return std::log(gaus3d(x,p) );
00056 }
00057
00058 double gaus10d(double * x, double *) {
00059 int i;
00060 double tmp = 0.;
00061 for (i=0; i<10; i++)
00062 tmp -= x[i] * x[i];
00063 return exp(tmp/2.);
00064 }
00065
00066 double gaus100d(double * x, double *) {
00067 int i;
00068 double tmp = 0.;
00069 for (i=0; i<100; i++)
00070 tmp -= x[i] * x[i];
00071 return exp(tmp/2.);
00072 }
00073
00074
00075 int testUnuran(TUnuran & unr, const std::string & method, const TUnuranMultiContDist & dist, TH3 * h1, const TH3 * href = 0, const int dim = 3 ) {
00076
00077
00078 assert (dim >=3);
00079
00080 bool ret = unr.Init(dist,method);
00081 if (!ret) {
00082 std::cerr << "Error initializing unuran with method " << unr.MethodName() << endl;
00083 return -1;
00084 }
00085
00086 h1->Reset();
00087
00088
00089
00090 TStopwatch w;
00091
00092 w.Start();
00093 std::vector<double> x(dim);
00094 for (int i = 0; i < n; ++i) {
00095 unr.SampleMulti( &x[0]);
00096 h1->Fill(x[0],x[1],x[2]);
00097 }
00098
00099 w.Stop();
00100 double time = w.CpuTime()*1.E9/n;
00101
00102 cout << "Time using Unuran " << unr.MethodName() << " \t=\t " << time << "\tns/call\t";
00103 if (href != 0) {
00104 double prob = href->Chi2Test(h1,"UU");
00105 double ksprob = href->KolmogorovTest(h1);
00106 std::cout << "\tChi2 Prob = "<< prob << "\tKS Prob = " << ksprob << std::endl;
00107
00108
00109 if (unr.MethodName() == "hitro") prob = ksprob;
00110 if (prob < 1.E-6 ) {
00111 std::cout << "\nChi2 Test failed ! " << std::endl;
00112 href->Chi2Test(h1,"UUP");
00113 return 1;
00114 }
00115 }
00116 else
00117 std::cout << std::endl;
00118
00119 return 0;
00120 }
00121
00122 int testGetRandom(TF3 * f, TH3 * h1, const TH3 * href = 0) {
00123
00124
00125
00126 h1->Reset();
00127
00128 TStopwatch w;
00129 w.Start();
00130 double x[3] = {0,0,0};
00131 for (int i = 0; i < n; ++i) {
00132 f->GetRandom3(x[0],x[1],x[2]);
00133 h1->Fill(x[0],x[1],x[2]);
00134 }
00135
00136 w.Stop();
00137 double time = w.CpuTime()*1.E9/n;
00138
00139
00140 if (href != 0) {
00141 double prob = href->Chi2Test(h1,"UU");
00142 double ksprob = href->KolmogorovTest(h1);
00143 std::cout << "Time using TF1::GetRandom() \t=\t " << time << "\tns/call \t\tChi2 Prob = "<< prob
00144 << "\tKS Prob = " << ksprob << std::endl;
00145 if (prob < 1E-06) {
00146 std::cout << "\tChi2 Test failed ! " << std::endl;
00147 href->Chi2Test(h1,"UUP");
00148 return 2;
00149 }
00150 }
00151 else
00152 std::cout << "Time using TF1::GetRandom() \t=\t " << time << "\tns/call\n";
00153 return 0;
00154 }
00155
00156
00157 int unuranMultiDim() {
00158
00159
00160 gErrorIgnoreLevel = 1001;
00161
00162
00163 if (useRandomSeed) gRandom->SetSeed(0);
00164
00165 gSystem->Load("libMathCore");
00166 gSystem->Load("libUnuran");
00167
00168
00169 n = 100000;
00170
00171
00172
00173 TH3D * h1 = new TH3D("h1","UNURAN gaussian 3D distribution",50,-10,10,50,-10,10,50,-10,10);
00174 TH3D * h2 = new TH3D("h2","TF1::GetRandom gaussian 3D distribution",50,-10,10,50,-10,10,50,-10,10);
00175
00176
00177 TH3D * h3 = new TH3D("h3","UNURAN truncated gaussian 3D distribution",50,-1,1,50,-1,1,50,-1,1);
00178 TH3D * h4 = new TH3D("h4","TF1::GetRandom truncated gaussian 3D distribution",50,-1,1,50,-1,1,50,-1,1);
00179
00180
00181 TF3 * f = new TF3("g3d",gaus3d,-10,10,-10,10,-10,10,4);
00182 double par[4] = {2,2,2,0.5};
00183 f->SetParameters(par);
00184 TF3 * flog = new TF3("logg3d",log_gaus3d,-10,10,-10,10,-10,10,4);
00185 flog->SetParameters(par);
00186
00187
00188
00189 std::cout << "Test using an undefined domain :\n\n";
00190
00191 std::cout << "Function Points used in GetRandom: [ " << f->GetNpx() << " , "
00192 << f->GetNpy() << " , " << f->GetNpz() << " ]" << std::endl;
00193 testGetRandom(f,h1);
00194
00195
00196
00197 int np = 100;
00198 f->SetNpx(np); f->SetNpy(np); f->SetNpz(np);
00199 std::cout << "Function Points used in GetRandom: [ " << f->GetNpx() << " , "
00200 << f->GetNpy() << " , " << f->GetNpz() << " ]" << std::endl;
00201
00202 testGetRandom(f,h2,h1);
00203
00204 *h1 = *h2;
00205 np = 200;
00206 f->SetNpx(np); f->SetNpy(np); f->SetNpz(np);
00207 std::cout << "Function Points used in GetRandom: [ " << f->GetNpx() << " , "
00208 << f->GetNpy() << " , " << f->GetNpz() << " ]" << std::endl;
00209
00210 testGetRandom(f,h2,h1);
00211
00212
00213
00214 TUnuranMultiContDist dist(f);
00215
00216
00217 TUnuran unr(gRandom,2);
00218
00219 int iret = 0;
00220 TH3 * href = new TH3D(*h2);
00221
00222
00223 std::string method = "vnrou";
00224 iret |= testUnuran(unr, method, dist, h1, href);
00225
00226
00227
00228 method = "hitro";
00229 iret |= testUnuran(unr, method, dist, h1, href);
00230
00231
00232 #ifdef USE_GIBBS
00233 method = "gibbs";
00234
00235 TUnuranMultiContDist logdist(flog,0,true);
00236 iret |= testUnuran(unr, method, logdist, h1, href);
00237 #endif
00238
00239
00240 cout << "\nTest setting the mode in Unuran distribution:\n\n";
00241 double m[3] = {0,0,0};
00242 dist.SetMode(m);
00243
00244 method = "vnrou";
00245 iret |= testUnuran(unr, method, dist, h1, href);
00246
00247 method = "hitro";
00248 iret |= testUnuran(unr, method, dist, h1, href);
00249
00250 #ifdef USE_GIBBS
00251 method = "gibbs";
00252 logdist.SetMode(m);
00253 iret |= testUnuran(unr, method, logdist, h1, href);
00254 #endif
00255
00256 cout << "\nTest truncated distribution:\n\n";
00257
00258 double xmin[3] = { -1, -1, -1 };
00259 double xmax[3] = { 1, 1, 1 };
00260
00261 f->SetRange(xmin[0],xmin[1],xmin[2],xmax[0],xmax[1],xmax[2]);
00262
00263 dist.SetDomain(xmin,xmax);
00264
00265 np = 30;
00266 f->SetNpx(np); f->SetNpy(np); f->SetNpz(np);
00267 std::cout << "Function Points used in GetRandom: [ " << f->GetNpx() << " , "
00268 << f->GetNpy() << " , " << f->GetNpz() << " ]" << std::endl;
00269 testGetRandom(f,h3);
00270
00271 href = h3;
00272 np = 100;
00273 f->SetNpx(np); f->SetNpy(np); f->SetNpz(np);
00274 std::cout << "Function Points used in GetRandom: [ " << f->GetNpx() << " , "
00275 << f->GetNpy() << " , " << f->GetNpz() << " ]" << std::endl;
00276 testGetRandom(f,h4,href);
00277 href = h4;
00278
00279 method = "vnrou";
00280 iret |= testUnuran(unr, method, dist, h3, href);
00281 method = "hitro";
00282 iret |= testUnuran(unr, method, dist, h3, href);
00283
00284
00285 #ifdef USE_GIBBS
00286 method = "gibbs";
00287 logdist.SetDomain(xmin,xmax);
00288 iret |= testUnuran(unr, method, logdist, h3, href);
00289 #endif
00290
00291
00292 TCanvas * c1 = new TCanvas("c1_unuranMulti","Multidimensional distribution",10,10,900,900);
00293 c1->Divide(2,2);
00294
00295 c1->cd(1); h1->Draw();
00296 c1->cd(2); h2->Draw();
00297 c1->cd(3); h3->Draw();
00298 c1->cd(4); h4->Draw();
00299
00300
00301 c1->Update();
00302
00303
00304
00305 TH3D * hrefN = new TH3D("hrefN","UNURAN gaussian ref N-Dim distribution (first 3 dim)",30,-3,3,30,-3,3,30,-3,3);
00306 TH3D * h10v = new TH3D("h10v","UNURAN gaussian N-Dim distribution (first 3 dim)",30,-3,3,30,-3,3,30,-3,3);
00307 TH3D * h10h = new TH3D("h10h","UNURAN gaussian N-Dim distribution (first 3 dim)",30,-3,3,30,-3,3,30,-3,3);
00308 TH3D * h100 = new TH3D("h100","UNURAN gaussian N-Dim distribution (first 3 dim)",30,-3,3,30,-3,3,30,-3,3);
00309
00310 int scale = 5;
00311
00312 double par2[4] = {1,1,1,0.};
00313 f->SetParameters(par2);
00314 TUnuranMultiContDist dist3(f);
00315 method = "vnrou";
00316 n/= scale;
00317 iret |= testUnuran(unr, method, dist3, hrefN);
00318
00319 cout << "\nTest 10 dimension: (be patient , it takes time....)\n\n";
00320
00321 TF1 * f10 = new TF1("g10d",gaus10d,-10,10,0);
00322 TUnuranMultiContDist dist10(f10,10);
00323
00324
00325 TCanvas * c2 = new TCanvas("c2_unuranMulti","Multidimensional distribution",100,10,900,900);
00326 c2->Divide(2,2);
00327
00328 c2->cd(1); hrefN->Draw();
00329
00330
00331 method = "vnrou";
00332 iret |= testUnuran(unr, method, dist10, h10v, hrefN,10);
00333 c2->cd(2); h10v->Draw();
00334 c2->Update();
00335
00336
00337 method = "hitro;thinning=20";
00338 iret |= testUnuran(unr, method, dist10, h10h, hrefN,10);
00339 c2->cd(3); h10h->Draw();
00340 c2->Update();
00341
00342
00343
00344 cout << "\nTest 100 dimension: ( be patient , it takes time....)\n\n";
00345 TF1 * f100 = new TF1("g100d",gaus100d,-10,10,0);
00346 TUnuranMultiContDist dist100(f100,100);
00347
00348
00349
00350 std::cout << "number of events to be generated = " << n << endl;
00351 method = "hitro;thinning=200";
00352 iret |= testUnuran(unr, method, dist100, h100, hrefN,100);
00353
00354
00355
00356
00357
00358 c2->cd(4); h100->Draw();
00359 c2->Update();
00360
00361
00362 if (iret != 0)
00363 std::cerr <<"\n\nUnuRan MultiDim Continous Distribution Test:\t Failed !!!!!!!\n" << std::endl;
00364 else
00365 std::cerr << "\n\nUnuRan MultiDim Continous Distribution Test:\t OK\n" << std::endl;
00366 return iret;
00367
00368 }
00369
00370 #ifndef __CINT__
00371 int main(int argc, char **argv)
00372 {
00373 int iret = 0;
00374 if (argc > 1) {
00375 TApplication theApp("App",&argc,argv);
00376 iret = unuranMultiDim();
00377 theApp.Run();
00378 }
00379 else
00380 iret = unuranMultiDim();
00381
00382 return iret;
00383 }
00384 #endif
00385