unuranHist.cxx

Go to the documentation of this file.
00001 //test  Unuran using as  empirical distribution. 
00002 // Use as input an histogram (using its buffer) or a TGraph and TGraph2D for multi-dimensional data
00003 //
00004 // run the test within ROOT (.x unuranHist.cxx+) or pass any extra parameter in the command line to get  
00005 // a graphics output  (./unuranHist 1) 
00006 // 
00007 
00008 #include "TUnuran.h"
00009 #include "TUnuranEmpDist.h"
00010 #include "TUnuranMultiContDist.h"
00011 
00012 #include "TH1.h"
00013 #include "TH2.h"
00014 #include "TH3.h"
00015 #include "TMath.h"
00016 #include "TF1.h"
00017 #include "TF2.h"
00018 #include "TF3.h"
00019 #include "TGraph.h"
00020 #include "TGraph2D.h"
00021 #include "TApplication.h"
00022 #include "TCanvas.h"
00023 #include "TRandom3.h"
00024 
00025 #include "TStopwatch.h"
00026 #include "TError.h"
00027 
00028 #include <iterator> 
00029 
00030 #include <iostream>
00031 
00032 
00033 int unuranHist() {
00034 
00035    // switch off printing of  info messages from chi2 test
00036    gErrorIgnoreLevel = 1001; 
00037 
00038 
00039    int nbin = 100;
00040    double xmin = -5; 
00041    double xmax = 20;
00042 
00043    int n = 100000; // number of events generated and of reference histo
00044    int ns = 10000;   // number of events from starting histo
00045 
00046    // h0 is reference histo
00047    TH1D * h0 = new TH1D("h0","Landau ref data",nbin,xmin,xmax);
00048    h0->FillRandom("landau",n); 
00049    
00050    // h1 is low statistic histo from where we generate numbers
00051    TH1D * h1 = new TH1D("h1","Landau source data",nbin,xmin,xmax);
00052    h1->SetBuffer(ns);
00053 
00054 
00055    h1->FillRandom("landau",ns); 
00056 
00057 //    const double * bf = h1->GetBuffer();
00058 //    std::ostream_iterator<double> oi(std::cout," ,  ");
00059 //    std::copy(bf, bf+ 2*h1->GetBufferLength() + 2, oi); 
00060 //    std::cout << std::endl << std::endl; 
00061  
00062    int iret = 0;
00063 
00064    std::cout << "Test Using 1D UnBinned data\n " << std::endl;
00065 
00066 
00067    TUnuran unr; 
00068    TUnuranEmpDist dist(h1);
00069 
00070    if (!unr.Init(dist)) return -1;
00071 
00072    TStopwatch w; 
00073    w.Start(); 
00074    TH1D * h2 = new TH1D("h2","Landau from Unuran unbin generation",nbin,xmin,xmax);
00075    for (int i = 0; i < n; ++i) { 
00076       h2->Fill( unr.Sample() );
00077    }
00078    w.Stop(); 
00079    double time = w.CpuTime()*1.E9/n; 
00080    std::cout << "Time using Unuran  " << unr.MethodName() << "   \t=\t " << time << "\tns/call" << std::endl;
00081 
00082    w.Start();
00083    TH1D * h3 = new TH1D("h3","h3",nbin,xmin,xmax);
00084    for (int i = 0; i < n; ++i) { 
00085       h3->Fill( h1->GetRandom() );
00086    }
00087    w.Stop(); 
00088    time = w.CpuTime()*1.E9/n;
00089    std::cout << "Time using TH1::GetRandom()  \t=\t " <<  time << "\tns/call" << std::endl;
00090 
00091 
00092    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << " (with h0)\t:\t";
00093    double prob = h2->Chi2Test(h1,"UUP");
00094    if (prob < 1.E-6 ) { 
00095       std::cerr << "Chi2 Test failed for UNURAN method " <<  unr.MethodName() << std::endl;
00096       iret = -1;
00097    }
00098    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << " (with ref)\t:\t";
00099    h2->Chi2Test(h0,"UUP");
00100 
00101    std::cout << "Test quality TH1::GetRandom (with h1) \t:\t"; 
00102    h3->Chi2Test(h1,"UUP");
00103    std::cout << "Test quality TH1::GetRandom (with ref) \t:\t"; 
00104    h3->Chi2Test(h0,"UUP");
00105    std::cout << "Comparison UnuRan-TH1::GetRandom \t:\t"; 
00106    h2->Chi2Test(h3,"UUP");
00107 
00108 
00109 
00110    TCanvas * c1 = new TCanvas("c1_unuranEmp","Empirical 1D distribution",10,10,800,800); 
00111    c1->Divide(1,2);
00112    c1->cd(1);
00113    h2->SetLineColor(kBlue);
00114    h2->Draw();
00115    h0->Draw("same");
00116    h3->SetLineColor(kRed);
00117    h3->Draw("same");
00118 
00119 
00120 
00121    // generate using binned data from h1
00122    std::cout << "\nTest Using Binned data\n " << std::endl;
00123 
00124    TUnuranEmpDist dist2(h1,false);
00125 
00126 
00127    if (!unr.Init(dist2) ) return -1;
00128 
00129 
00130 
00131    w.Start(); 
00132    TH1D * h4 = new TH1D("h4","Landau from Unuran binned generation",nbin,xmin,xmax);
00133    for (int i = 0; i < n; ++i) { 
00134       h4->Fill( unr.Sample() );
00135    }
00136    w.Stop(); 
00137    time = w.CpuTime()*1.E9/n; 
00138    std::cout << "Time using Unuran  " << unr.MethodName() << "   \t=\t " << time << "\tns/call" << std::endl;
00139 
00140    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << " (with h1)\t:\t"; 
00141    double prob2 = h4->Chi2Test(h1,"UUP");
00142    if (prob2 < 1.E-6) { 
00143       std::cerr << "Chi2 Test failed for UNURAN method " <<  unr.MethodName() << std::endl;
00144       iret = -2;
00145    }
00146    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << " (with ref)\t:\t"; 
00147    h4->Chi2Test(h0,"UUP");
00148 
00149    c1->cd(2);
00150    h4->SetLineColor(kBlue);
00151    h4->Draw();
00152    h0->Draw("same");
00153 
00154    return iret; 
00155 }
00156 
00157 double gaus2d(double *x, double *p) { 
00158 
00159    double sigma_x = p[0]; 
00160    double sigma_y = p[1];
00161    double rho = p[2]; 
00162    double u = x[0] / sigma_x ;
00163    double v = x[1] / sigma_y ;
00164    double c = 1 - rho*rho ;
00165    double result = (1 / (2 * TMath::Pi() * sigma_x * sigma_y * sqrt(c))) 
00166       * exp (-(u * u - 2 * rho * u * v + v * v ) / (2 * c));
00167    return result;
00168 }
00169 
00170 double gaus3d(double *x, double *p) { 
00171 
00172    double sigma_x = p[0]; 
00173    double sigma_y = p[1];
00174    double sigma_z = p[2];
00175    double rho = p[3]; 
00176    double u = x[0] / sigma_x ;
00177    double v = x[1] / sigma_y ;
00178    double w = x[2] / sigma_z ;
00179    double c = 1 - rho*rho ;
00180    double result = (1 / (2 * TMath::Pi() * sigma_x * sigma_y * sigma_z * sqrt(c))) 
00181       * exp (-(u * u - 2 * rho * u * v + v * v + w*w) / (2 * c));
00182    return result;
00183 }
00184 
00185 
00186 int unuranGraf() { 
00187 
00188    int iret = 0; 
00189 
00190    std::cout << "\nTest Using 2D UnBinned data\n " << std::endl;
00191 
00192    // generate graf with x-y data
00193    const int Ndata = 100000;
00194    const int n = 100000;
00195 
00196    TF2 * f2 = new TF2("g2d",gaus2d,-10,10,-10,10,3);
00197    double p2[3] = {1,1,0.5};
00198    f2->SetParameters(p2);
00199 
00200    TUnuran unr(gRandom,2); 
00201    if (!unr.Init(TUnuranMultiContDist(f2),"vnrou") ) return -1; 
00202    // generate 2d data
00203    double xx[2];
00204    // generate ref data
00205    TH2D * href = new TH2D("href2d","UNURAN reference 2D gauss data",100,-5,5,100,-5,5);  
00206    for (int i = 0; i < n; ++i) { 
00207       unr.SampleMulti(xx);
00208       href->Fill(xx[0],xx[1]);
00209    }
00210 
00211    // create a graph with source xy data
00212    TGraph * gr = new TGraph();
00213 
00214    for (int i = 0; i< Ndata; ++i) { 
00215       unr.SampleMulti( xx); 
00216       gr->SetPoint(i,xx[0],xx[1]);
00217    }
00218 
00219    TH2D * h2 = new TH2D("h2d","UNURAN generated 2D gauss data",100,-5,5,100,-5,5);  
00220 
00221    TH1D * hx = new TH1D("hx","x gen variable",100,-5,5);  
00222    TH1D * hy = new TH1D("hy","y gen variable",100,-5,5);  
00223    
00224    // now generate random points from the graph
00225    TUnuranEmpDist dist(Ndata,gr->GetX(), gr->GetY() ); 
00226    if (!unr.Init(dist) ) { 
00227       std::cerr << "Initialization failed for method " << unr.MethodName() << std::endl;
00228       return -1;
00229    }
00230 
00231    
00232    TStopwatch w; 
00233    w.Start(); 
00234    for (int i = 0; i < n; ++i) { 
00235       unr.SampleMulti(xx);
00236       h2->Fill(xx[0],xx[1]);
00237       hx->Fill(xx[0]);
00238       hy->Fill(xx[1]);
00239    }
00240    w.Stop(); 
00241    double time = w.CpuTime()*1.E9/n; 
00242    std::cout << "Time using Unuran  " << unr.MethodName() << "   \t=\t " << time << "\tns/call" << std::endl;
00243 
00244    TCanvas * c1 = new TCanvas("c1_unuranEmp2D","Empirical Multidim distribution",300,10,800,400); 
00245 
00246    c1->Divide(2,1);
00247    c1->cd(1);
00248    gr->Draw("AP");
00249    c1->cd(2);
00250    h2->Draw("surf"); 
00251    //f2->Draw("same");
00252 
00253    TCanvas * c2 = new TCanvas("c2d","Empirical Multidim distribution",300,100,800,400); 
00254    c2->Divide(2,1);
00255    c2->cd(1);
00256    hx->Draw();
00257    c2->cd(2);
00258    hy->Draw(); 
00259 
00260    // apply chi2 test to href
00261    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << "\t\t:\t"; 
00262    double prob = href->Chi2Test(h2,"UUP");
00263    if (prob < 1.E-6) { 
00264       std::cerr << "Chi2 Test failed for UNURAN method " <<  unr.MethodName() << std::endl;
00265       iret = -2;
00266    }
00267 
00268    return iret; 
00269 }
00270 
00271 int unuranGraf2D() { 
00272    int iret = 0; 
00273 
00274    std::cout << "\nTest Using 3D UnBinned data\n " << std::endl;
00275 
00276    // generate graf with x-y data
00277    const int Ndata = 100000;
00278    const int n = 100000;
00279 
00280    TF3 * f3 = new TF3("g3d",gaus3d,-10,10,-10,10,-10,10,4);
00281    double p[4] = {1,1,1,0.5};
00282    f3->SetParameters(p);
00283 
00284    TUnuran unr(gRandom,2); 
00285    if (!unr.Init(TUnuranMultiContDist(f3),"vnrou") ) return -1; 
00286    // generate 3d data
00287    double xx[3];
00288    // generate ref data
00289    TH3D * href = new TH3D("href3d","UNURAN reference 3D gauss data",50,-5,5,50,-5,5,50,-5,5);  
00290 
00291    for (int i = 0; i < n; ++i) { 
00292       unr.SampleMulti(xx);
00293       href->Fill(xx[0],xx[1],xx[2]);
00294    }
00295 
00296    // create a graph with xy data
00297    TGraph2D * gr = new TGraph2D();
00298    for (int i = 0; i< Ndata; ++i) { 
00299       unr.SampleMulti( xx); 
00300       gr->SetPoint(i,xx[0],xx[1],xx[2]);
00301    }
00302 
00303    TH3D * h3 = new TH3D("h3d","UNURAN generated 3D gauss data",50,-5,5,50,-5,5,50,-5,5);  
00304 
00305    TH1D * hx = new TH1D("hx3","x gen variable",100,-5,5);  
00306    TH1D * hy = new TH1D("hy3","y gen variable",100,-5,5);  
00307    TH1D * hz = new TH1D("hz3","z gen variable",100,-5,5);  
00308 
00309    // now generate random points from the graph
00310    TUnuranEmpDist dist(Ndata,gr->GetX(), gr->GetY(), gr->GetZ() ); 
00311    if (!unr.Init(dist) ) { 
00312       std::cerr << "Initialization failed for method " << unr.MethodName() << std::endl;
00313       return -1;
00314    } 
00315 
00316    
00317    TStopwatch w; 
00318    w.Start(); 
00319    for (int i = 0; i < n; ++i) { 
00320       unr.SampleMulti(xx);
00321       h3->Fill(xx[0],xx[1],xx[2]);
00322       hx->Fill(xx[0]);
00323       hy->Fill(xx[1]);
00324       hz->Fill(xx[2]);
00325    }
00326    w.Stop(); 
00327    double time = w.CpuTime()*1.E9/n; 
00328    std::cout << "Time using Unuran  " << unr.MethodName() << "   \t=\t " << time << "\tns/call" << std::endl;
00329 
00330    TCanvas * c3 = new TCanvas("c3d","Empirical Multidim distribution",300,600,800,400); 
00331    c3->Divide(3,1);
00332    c3->cd(1);
00333    hx->Draw();
00334    c3->cd(2);
00335    hy->Draw(); 
00336    c3->cd(3);
00337    hz->Draw(); 
00338 
00339    TCanvas * c1 = new TCanvas("c1_unuranEmp3D","Empirical Multidim distribution",400,500,900,400); 
00340    c1->Divide(2,1); 
00341    c1->cd(1);
00342    gr->Draw("AP");
00343    c1->cd(2);
00344    h3->Draw("surf1"); 
00345    f3->Draw("same");
00346 
00347    // apply chi2 test to href
00348    std::cout << "\nTest quality UNURAN  " << unr.MethodName() << "\t\t:\t"; 
00349    double prob = href->Chi2Test(h3,"UUP");
00350    if (prob < 1.E-6) { 
00351       std::cerr << "Chi2 Test failed for UNURAN method " <<  unr.MethodName() << std::endl;
00352       iret = -2;
00353    }
00354 
00355    return iret; 
00356 }
00357 
00358 #ifndef __CINT__
00359 int main(int argc, char **argv)
00360 {
00361    int iret = 0; 
00362    if (argc > 1) { 
00363       TApplication theApp("App",&argc,argv);
00364       iret |=  unuranHist();
00365       iret |=  unuranGraf();
00366       iret |=  unuranGraf2D();
00367       theApp.Run();
00368    } 
00369    else {
00370       iret |=  unuranHist();
00371       iret |=  unuranGraf();
00372       iret |=  unuranGraf2D();
00373    }
00374    
00375    if (iret != 0) 
00376       std::cerr <<"\n\nUnuRan Empirical Distribution Test:\t  Failed !!!!!!!\n" << std::endl;
00377    else 
00378       std::cout << "\n\nUnuRan  Empirical Distribution Test:\t OK\n" << std::endl;
00379    return iret; 
00380 }
00381 #endif

Generated on Tue Jul 5 14:43:43 2011 for ROOT_528-00b_version by  doxygen 1.5.1