![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
Random number generators are descendants of rnga . While the base object rnga is created to allow user-defined random number generators, the only random number generator presently included are from GSL. Because speed is frequently an issue, some design choices were made to ensure fast execution
This examples creates a Markov chain for an aritrary distribution using the Metropolis-Hastings algorithm. It generates the one-dimensional distribution
for .
/* Example: ex_markov.cpp ------------------------------------------------------------------- An example to demonstrate generation of an arbitrary distribution through the creation of a Markov chain using the Metropolis algorithm. */ #include <iostream> #include <cmath> #include <o2scl/test_mgr.h> #include <o2scl/hist.h> #include <o2scl/gsl_rnga.h> #include <o2scl/expect_val.h> #include <o2scl/gsl_inte_qag.h> using namespace std; using namespace o2scl; /* An unusual two-peaked probability distribution. This distribution is is strongly bimodal, and thus the algorithm has to step over a region where the probability distribution is small. */ double function(double x) { return exp(-x*x)*(sin(x-1.4)+1.0); } // Same function but now weighted with x double function2(double x) { return x*exp(-x*x)*(sin(x-1.4)+1.0); } int main(int argc, char *argv[]) { test_mgr t; t.set_output_level(1); cout.setf(ios::scientific); // The square root of the number of iterations static const size_t N=1000; // The stepsize for the Metropolis algorithm. If this stepsize was // too small, the algorithm would fail to step between the two // peaks and give incorrect results. double step=1.0; // The random number generator gsl_rnga gr; if (argc>=2) { gr.set_seed(stoi(argv[1])); } // Construct a histogrmam of the distribution. Initialize the grid // spacing with a unform grid between -5 and 5. hist h; h.set_bins(uniform_grid_end_width<>(-5.0,5.0,0.1)); // Compute the average value of the distribution, and use // a scalar_ev object to estimate the uncertainty static const size_t M=20; scalar_ev sev(M,N*N/M); // An initial point and the initial weight double x_old=0.01; double w_old=function(x_old); for(size_t k=0;k<N;k++) { for(size_t i=0;i<N;i++) { // Perform a step and compute the new weight double r=gr.random(); double x_new=x_old+r*step*2.0-step; double w_new=function(x_new); // The Metrpolis algorithm double r2=gr.random(); if (r2<w_new/w_old && x_new>-5.0 && x_new<5.0) { x_old=x_new; w_old=w_new; } // Update the histogram h.update(x_old); // Update the average value sev.add(x_old); } } // Normalize by the value of the distribution at x=1.05 double norm=h.get_wgt(1.05)/function(1.05); // Output the generated distribution and compare to the real // distribution. for(size_t i=0;i<h.size();i++) { double xx=h.get_rep_i(i); // Limit output to every third point if (i%3==0) { cout.width(2); cout << i << " "; cout.setf(ios::showpos); cout << xx << " "; cout.unsetf(ios::showpos); cout << h.get_wgt_i(i)/norm << " " << function(xx) << endl; } // Only test regions where we'll have enough statistics if (function(xx)>5.0e-3) { t.test_rel(h.get_wgt_i(i)/norm,function(xx),1.0e-1,"dist"); } } cout << endl; // Output the average value and its uncertainty double avg, std_dev, avg_err; sev.current_avg(avg,std_dev,avg_err); cout << avg << " " << avg_err << endl; // Compare to result from direct integration funct_fptr ff(function); funct_fptr ff2(function2); gsl_inte_qag<> qag; double exact=qag.integ(ff2,-5.0,5.0)/qag.integ(ff,-5.0,5.0); cout << exact << endl; t.test_abs(avg,exact,3.0*avg_err,"Average value"); t.report(); return 0; } // End of example
Typical output is
0 -4.950000e+00 0.000000e+00 2.131525e-11 3 -4.650000e+00 0.000000e+00 5.009023e-10 6 -4.350000e+00 0.000000e+00 9.131547e-09 9 -4.050000e+00 0.000000e+00 1.309343e-07 12 -3.750000e+00 4.095279e-06 1.488688e-06 15 -3.450000e+00 1.228584e-05 1.348287e-05 18 -3.150000e+00 5.323863e-05 9.747130e-05 21 -2.850000e+00 5.487674e-04 5.624059e-04 24 -2.550000e+00 2.706979e-03 2.584240e-03 27 -2.250000e+00 9.587048e-03 9.410936e-03 30 -1.950000e+00 2.736465e-02 2.693191e-02 33 -1.650000e+00 5.785401e-02 5.970011e-02 36 -1.350000e+00 1.004859e-01 9.993669e-02 39 -1.050000e+00 1.204053e-01 1.202766e-01 42 -7.500000e-01 9.158682e-02 9.293227e-02 45 -4.500000e-01 3.099307e-02 3.162602e-02 48 -1.500000e-01 6.757210e-04 2.114248e-04 51 +1.500000e-01 4.999517e-02 4.988035e-02 54 +4.500000e-01 1.519635e-01 1.523810e-01 57 +7.500000e-01 2.238029e-01 2.249580e-01 60 +1.050000e+00 2.181842e-01 2.181842e-01 63 +1.350000e+00 1.529259e-01 1.535435e-01 66 +1.650000e+00 8.190148e-02 8.196725e-02 69 +1.950000e+00 3.340519e-02 3.397864e-02 72 +2.250000e+00 1.076239e-02 1.108511e-02 75 +2.550000e+00 3.153365e-03 2.868544e-03 78 +2.850000e+00 7.576266e-04 5.914089e-04 81 +3.150000e+00 5.733391e-05 9.733109e-05 84 +3.450000e+00 4.095279e-06 1.278395e-05 87 +3.750000e+00 0.000000e+00 1.336916e-06 90 +4.050000e+00 0.000000e+00 1.107648e-07 93 +4.350000e+00 0.000000e+00 7.207155e-09 96 +4.650000e+00 0.000000e+00 3.628586e-10 99 +4.950000e+00 0.000000e+00 1.376924e-11 46 tests performed. All tests passed.
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).