unuranDemo.C

Go to the documentation of this file.
00001 // Example macro to show unuran capabilities 
00002 // The results are compared with what is obtained using TRandom or TF1::GetRandom
00003 // The macro is divided in 3 parts: 
00004 //    - testStringAPI         :  show how to use string API of UNURAN to generate Gaussian random numbers
00005 //    - testDistr1D           :  show how to pass a 1D distribution object to UNURAN to generate numbers 
00006 //                               according to the given distribution object
00007 //    - testDistrMultiDIm     :  show how to pass a multidimensional distribution object to UNURAN 
00008 //                               
00009 //
00010 // To execute the macro type in: 
00011 //
00012 // root[0]: gSystem->Load("libMathCore");
00013 // root[0]: gSystem->Load("libUnuran");
00014 // root[0]: .x  unuranDemo.C+
00015 //
00016 //Author: Lorenzo Moneta
00017 
00018 
00019 #include "TStopwatch.h"
00020 
00021 #include "TUnuran.h"
00022 #include "TUnuranContDist.h"
00023 #include "TUnuranMultiContDist.h"
00024 #include "TUnuranDiscrDist.h"
00025 #include "TUnuranEmpDist.h"
00026 
00027 #include "TH1.h"
00028 #include "TH3.h"
00029 #include "TF3.h"
00030 #include "TMath.h"
00031 
00032 #include "TRandom2.h"
00033 #include "TSystem.h"
00034 #include "TStyle.h"
00035 
00036 #include "TApplication.h"
00037 #include "TCanvas.h"
00038 
00039 #ifndef __CINT__  // need to exclude to avoid CINT re-defining them
00040 #include "Math/ProbFunc.h"
00041 #include "Math/DistFunc.h"
00042 #endif
00043 
00044 #include <iostream> 
00045 #include <cassert>
00046 
00047 using std::cout; 
00048 using std::endl; 
00049 
00050 // number of distribution generated points
00051 #define NGEN 1000000  
00052 
00053 int izone = 0; 
00054 TCanvas * c1 = 0; 
00055 
00056 // test using UNURAN string interface
00057 void testStringAPI() {
00058 
00059    TH1D * h1 = new TH1D("h1G","gaussian distribution from Unuran",100,-10,10);
00060    TH1D * h2 = new TH1D("h2G","gaussian distribution from TRandom",100,-10,10);
00061 
00062    cout << "\nTest using UNURAN string API \n\n";
00063 
00064 
00065    TUnuran unr; 
00066    if (! unr.Init( "normal()", "method=arou") ) {
00067       cout << "Error initializing unuran" << endl;
00068       return;
00069    }
00070 
00071    int n = NGEN;
00072    TStopwatch w; 
00073    w.Start(); 
00074 
00075    for (int i = 0; i < n; ++i) {
00076        double x = unr.Sample(); 
00077        h1->Fill(  x ); 
00078    }
00079 
00080    w.Stop(); 
00081    cout << "Time using Unuran method " << unr.MethodName() << "\t=\t " << w.CpuTime() << endl;
00082 
00083 
00084    // use TRandom::Gaus
00085    w.Start();
00086    for (int i = 0; i < n; ++i) {
00087       double x = gRandom->Gaus(0,1); 
00088        h2->Fill(  x ); 
00089    }
00090 
00091    w.Stop(); 
00092    cout << "Time using TRandom::Gaus  \t=\t " << w.CpuTime() << endl;
00093 
00094    assert(c1 != 0);
00095    c1->cd(++izone);
00096    h1->Draw();
00097    c1->cd(++izone);
00098    h2->Draw();
00099 
00100 }
00101 
00102 
00103 
00104 double distr(double *x, double *p) { 
00105    return ROOT::Math::breitwigner_pdf(x[0],p[0],p[1]); 
00106 }
00107 
00108 double cdf(double *x, double *p) { 
00109    return ROOT::Math::breitwigner_cdf(x[0],p[0],p[1]); 
00110 }
00111  
00112 // test of unuran passing as input a distribution object( a BreitWigner) distribution 
00113 void testDistr1D() { 
00114 
00115    cout << "\nTest 1D Continous distributions\n\n";
00116 
00117    TH1D * h1 = new TH1D("h1BW","Breit-Wigner distribution from Unuran",100,-10,10);
00118    TH1D * h2 = new TH1D("h2BW","Breit-Wigner distribution from GetRandom",100,-10,10);
00119 
00120    
00121 
00122    TF1 * f = new TF1("distrFunc",distr,-10,10,2); 
00123    double par[2] = {1,0};  // values are gamma and mean 
00124    f->SetParameters(par); 
00125 
00126    TF1 * fc = new TF1("cdfFunc",cdf,-10,10,2); 
00127    fc->SetParameters(par); 
00128 
00129    // create Unuran 1D distribution object 
00130    TUnuranContDist dist(f); 
00131    // to use a different random number engine do: 
00132    TRandom2 * random = new TRandom2();
00133    int logLevel = 2;
00134    TUnuran unr(random,logLevel); 
00135 
00136    // select unuran method for generating the random numbers
00137    std::string method = "tdr";
00138    //std::string method = "method=auto";
00139    // "method=hinv" 
00140    // set the cdf for some methods like hinv that requires it
00141    // dist.SetCdf(fc); 
00142 
00143    //cout << "unuran method is " << method << endl;
00144 
00145    if (!unr.Init(dist,method) ) { 
00146       cout << "Error initializing unuran" << endl;
00147       return;
00148    } 
00149 
00150 
00151 
00152    TStopwatch w; 
00153    w.Start(); 
00154 
00155    int n = NGEN;
00156    for (int i = 0; i < n; ++i) {
00157       double x = unr.Sample(); 
00158       h1->Fill(  x ); 
00159    }
00160 
00161    w.Stop(); 
00162    cout << "Time using Unuran method " << unr.MethodName() << "\t=\t " << w.CpuTime() << endl;
00163 
00164    w.Start();
00165    for (int i = 0; i < n; ++i) {
00166       double x = f->GetRandom(); 
00167       h2->Fill(  x ); 
00168    }
00169 
00170    w.Stop(); 
00171    cout << "Time using TF1::GetRandom()  \t=\t " << w.CpuTime() << endl;
00172 
00173    c1->cd(++izone);
00174    h1->Draw();
00175 
00176    c1->cd(++izone);
00177    h2->Draw();
00178    
00179    std::cout << " chi2 test of UNURAN vs GetRandom generated histograms:  " << std::endl;
00180    h1->Chi2Test(h2,"UUP");
00181 
00182 }
00183 
00184 // 3D gaussian distribution
00185 double gaus3d(double *x, double *p) { 
00186 
00187    double sigma_x = p[0]; 
00188    double sigma_y = p[1];
00189    double sigma_z = p[2];
00190    double rho = p[2]; 
00191    double u = x[0] / sigma_x ;
00192    double v = x[1] / sigma_y ;
00193    double w = x[2] / sigma_z ;
00194    double c = 1 - rho*rho ;
00195    double result = (1 / (2 * TMath::Pi() * sigma_x * sigma_y * sigma_z * sqrt(c))) 
00196       * exp (-(u * u - 2 * rho * u * v + v * v + w*w) / (2 * c));
00197    return result;
00198 }
00199 
00200 // test of unuran passing as input a mluti-dimension distribution object 
00201 void testDistrMultiDim() { 
00202    
00203    cout << "\nTest Multidimensional distributions\n\n";
00204 
00205    TH3D * h1 = new TH3D("h13D","gaussian 3D distribution from Unuran",50,-10,10,50,-10,10,50,-10,10);
00206    TH3D * h2 = new TH3D("h23D","gaussian 3D distribution from GetRandom",50,-10,10,50,-10,10,50,-10,10);
00207 
00208 
00209 
00210    TF3 * f = new TF3("g3d",gaus3d,-10,10,-10,10,-10,10,3); 
00211    double par[3] = {2,2,0.5}; 
00212    f->SetParameters(par); 
00213 
00214 
00215 
00216    TUnuranMultiContDist dist(f); 
00217    TUnuran unr(gRandom);
00218    //std::string method = "method=vnrou";
00219 
00220    //std::string method = "method=hitro;use_boundingrectangle=false "; 
00221    std::string method = "hitro"; 
00222    if ( !  unr.Init(dist,method) ) { 
00223       cout << "Error initializing unuran" << endl;
00224       return;
00225    } 
00226 
00227    TStopwatch w; 
00228    w.Start(); 
00229 
00230    double x[3]; 
00231    for (int i = 0; i < NGEN; ++i) {  
00232       unr.SampleMulti(x);
00233       h1->Fill(x[0],x[1],x[2]);
00234    }
00235 
00236    w.Stop(); 
00237    cout << "Time using Unuran method " << unr.MethodName() << "\t=\t\t " << w.CpuTime() << endl;
00238 
00239    assert(c1 != 0);
00240    c1->cd(++izone); 
00241    h1->Draw();
00242 
00243 
00244    // need to set a reasanable number of points in TF1 to get accettable quality from GetRandom to 
00245    int np = 200;
00246    f->SetNpx(np);
00247    f->SetNpy(np);
00248    f->SetNpz(np);
00249 
00250    w.Start();
00251    for (int i = 0; i < NGEN; ++i) { 
00252       f->GetRandom3(x[0],x[1],x[2]);
00253       h2->Fill(x[0],x[1],x[2]); 
00254    }
00255 
00256    w.Stop(); 
00257    cout << "Time using TF1::GetRandom  \t\t=\t " << w.CpuTime() << endl;
00258 
00259    
00260    c1->cd(++izone);
00261    h2->Draw();
00262 
00263    std::cout << " chi2 test of UNURAN vs GetRandom generated histograms:  " << std::endl;
00264    h1->Chi2Test(h2,"UUP");
00265 
00266 }
00267 //_____________________________________________
00268 //
00269 // example of discrete distributions 
00270 
00271 double poisson(double * x, double * p) { 
00272    return ROOT::Math::poisson_pdf(int(x[0]),p[0]);
00273 }
00274 
00275 void testDiscDistr() { 
00276 
00277    cout << "\nTest Discrete distributions\n\n";
00278 
00279    TH1D * h1 = new TH1D("h1PS","Unuran Poisson prob",20,0,20);
00280    TH1D * h2 = new TH1D("h2PS","Poisson dist from TRandom",20,0,20);
00281 
00282    double mu = 5; 
00283 
00284    TF1 * f = new TF1("fps",poisson,1,0,1);
00285    f->SetParameter(0,mu);
00286 
00287    TUnuranDiscrDist dist2 = TUnuranDiscrDist(f);
00288    TUnuran unr;
00289 
00290    // dari method (needs also the mode and pmf sum)
00291    dist2.SetMode(int(mu) );
00292    dist2.SetProbSum(1.0);
00293    bool ret = unr.Init(dist2,"dari");
00294    if (!ret) return;
00295 
00296    TStopwatch w; 
00297    w.Start(); 
00298 
00299    int n = NGEN;
00300    for (int i = 0; i < n; ++i) {
00301       int k = unr.SampleDiscr(); 
00302       h1->Fill( double(k) ); 
00303    }
00304 
00305    w.Stop(); 
00306    cout << "Time using Unuran method " << unr.MethodName() << "\t=\t\t " << w.CpuTime() << endl;
00307 
00308    w.Start();
00309    for (int i = 0; i < n; ++i) {
00310       h2->Fill(  gRandom->Poisson(mu) ); 
00311    }
00312    cout << "Time using TRandom::Poisson " << "\t=\t\t " << w.CpuTime() << endl;
00313 
00314    c1->cd(++izone);
00315    h1->SetMarkerStyle(20);
00316    h1->Draw("E");
00317    h2->Draw("same");
00318    
00319    std::cout << " chi2 test of UNURAN vs TRandom generated histograms:  " << std::endl;
00320    h1->Chi2Test(h2,"UUP");
00321 
00322 }
00323  
00324 //_____________________________________________
00325 //
00326 // example of empirical distributions 
00327    
00328 void testEmpDistr() { 
00329 
00330 
00331    cout << "\nTest Empirical distributions using smoothing\n\n";
00332 
00333    // start with a set of data 
00334    // for example 1000 two-gaussian data
00335    const int Ndata = 1000;
00336    double x[Ndata];
00337    for (int i = 0; i < Ndata; ++i) {
00338       if (i < 0.5*Ndata ) 
00339          x[i] = gRandom->Gaus(-1.,1.);
00340       else 
00341          x[i] = gRandom->Gaus(1.,3.);
00342    }
00343 
00344    TH1D * h0 = new TH1D("h0Ref","Starting data",100,-10,10);
00345    TH1D * h1 = new TH1D("h1Unr","Unuran unbin Generated data",100,-10,10);
00346    TH1D * h1b = new TH1D("h1bUnr","Unuran bin Generated data",100,-10,10);
00347    TH1D * h2 = new TH1D("h2GR","Data from TH1::GetRandom",100,-10,10);
00348 
00349    h0->FillN(Ndata,x,0,1); // fill histogram with starting data
00350       
00351    TUnuran unr; 
00352    TUnuranEmpDist dist(x,x+Ndata,1);
00353 
00354 
00355    TStopwatch w; 
00356    int n = NGEN;
00357 
00358    w.Start(); 
00359    if (!unr.Init(dist)) return;
00360    for (int i = 0; i < n; ++i) {
00361       h1->Fill( unr.Sample() ); 
00362    }
00363 
00364    w.Stop(); 
00365    cout << "Time using Unuran unbin  " << unr.MethodName() << "\t=\t\t " << w.CpuTime() << endl;
00366 
00367    TUnuranEmpDist binDist(h0);
00368 
00369    w.Start(); 
00370    if (!unr.Init(binDist)) return;
00371    for (int i = 0; i < n; ++i) {
00372       h1b->Fill( unr.Sample() ); 
00373    }
00374    w.Stop(); 
00375    cout << "Time using Unuran bin  " << unr.MethodName() << "\t=\t\t " << w.CpuTime() << endl;
00376 
00377    w.Start();
00378    for (int i = 0; i < n; ++i) {
00379       h2->Fill(  h0->GetRandom() ); 
00380    }
00381    cout << "Time using TH1::GetRandom " << "\t=\t\t " << w.CpuTime() << endl;
00382 
00383    c1->cd(++izone);
00384 
00385    h2->Draw();
00386    h1->SetLineColor(kRed);
00387    h1->Draw("same");
00388    h1b->SetLineColor(kBlue);
00389    h1b->Draw("same");
00390 
00391    
00392 }
00393 
00394 
00395 
00396 void unuranDemo() { 
00397 
00398    //gRandom->SetSeed(0);
00399 
00400    // load libraries
00401    gSystem->Load("libMathCore");
00402    gSystem->Load("libUnuran");
00403 
00404    // create canvas 
00405 
00406    c1 = new TCanvas("c1_unuranMulti","Multidimensional distribution",10,10,1000,1000); 
00407    c1->Divide(2,4);
00408    gStyle->SetOptFit();
00409 
00410    testStringAPI(); 
00411    c1->Update();
00412    testDistr1D();
00413    c1->Update();
00414    testDistrMultiDim();
00415    c1->Update();
00416    testDiscDistr();
00417    c1->Update();
00418    testEmpDistr();
00419    c1->Update();
00420 
00421 
00422 }

Generated on Tue Jul 5 15:45:12 2011 for ROOT_528-00b_version by  doxygen 1.5.1