All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Monte Carlo Integration

Monte Carlo integration is performed by descendants of mcarlo_inte (o2scl::mcarlo_plain, o2scl::mcarlo_miser, and o2scl::mcarlo_vegas). These routines are generally superior to the direct methods for integrals over regions with large numbers of spatial dimensions.

Monte Carlo integration example

This example computes the integral

\[ \int_0^{\pi} \int_0^{\pi} \int_0^{\pi} \frac{1}{\pi^3} \left(1-\cos x \cos y \cos z\right)^{-1} \]

and compares the result with the exact value, 1.3932039296.

/* Example: ex_mcarlo.cpp
-------------------------------------------------------------------
An example to demonstrate multidimensional integration
*/
#include <boost/numeric/ublas/vector.hpp>
#include <o2scl/test_mgr.h>
#include <o2scl/multi_funct.h>
#include <o2scl/inte_multi_comp.h>
#include <o2scl/inte_qng_gsl.h>
#include <o2scl/mcarlo_vegas.h>
/// For M_PI
#include <gsl/gsl_math.h>
using namespace std;
using namespace o2scl;
double test_fun(size_t nv, const ubvector &x) {
double y=1.0/(1.0-cos(x[0])*cos(x[1])*cos(x[2]))/M_PI/M_PI/M_PI;
return y;
}
int main(void) {
cout.setf(ios::scientific);
double exact=1.3932039296;
double res;
double err;
ubvector a(3), b(3);
a[0]=0.0; b[0]=M_PI;
a[1]=0.0; b[1]=M_PI;
a[2]=0.0; b[2]=M_PI;
multi_funct11 tf=test_fun;
gm.n_points=100000;
gm.minteg_err(tf,3,a,b,res,err);
cout << res << " " << exact << " " << (res-exact)/err << endl;
t.test_rel(res,exact,err*10.0,"O2scl");
t.report();
return 0;
}
// End of example

Analysis of results from numerical simulations

The base o2scl::expval_base and its children form a set of classes useful for recording the outputs of several iterations of a numerical simulation, and summarizing the average, standard deviation, and error in the average over all iterations. After each iteration, some measurement is performed, and that measurement is added to the class with an add() functions (e.g. o2scl::expval_scalar::add() ).

Autocorrelations are common in numerical simulations, and these classes allow the user to assess the impact of correlations by breaking the measurements up into blocks.

Blocks are filled one by one, moving to the next block when the requested of measurements (n_per_block) in the previous block has been provided. If more measurements are given than the number of blocks times n_per_block, then the number of filled blocks is halved, each block is filled with twice as much data, and n_per_block is doubled. The blocks are always combined in such a way as to preserve the original ordering, i.e. the first block is combined with the second, the third with the fourth, and so on. Subsequent measurements are added to newly empty blocks.

The constructor for these classes need as input the number of blocks to record (at least one) and the number of measurements per block (also at least one). If either of these is specified to be zero then the error handler is called.

In order to properly handle autocorrelations, The number of measurements over which autocorrelation occurs should be less than the number of measurements per block.

Current averages are always reported as long as at least one measurement has been added, though not all data is always used if there are incomplete blocks or the blocks have not yet been filled evenly. Two functions current_avg() and current_avg_stats() are provided in children. The latter provides the number of blocks and number of measurements per block used in the currently reported statistics. If only one (full or partially full) block is available, the standard deviation and uncertainty on the average are reported as zero. If no data is available, then calls to current_avg() will call the error handler. See, e.g. o2scl::expval_scalar::current_avg() and o2scl::expval_scalar::current_avg_stats().

Many of the o2scl::expval_base children also support HDF I/O.

Note that the traditional Monte Carlo integration routines (o2scl::mcarlo_plain, o2scl::mcarlo_miser, and o2scl::mcarlo_vegas) have no significant autocorrelations by construction, so long as the random number generator does not have any.

Generate an arbitrary distribution

This examples creates a Markov chain for an aritrary distribution using the Metropolis-Hastings algorithm. It generates the one-dimensional distribution

\[ f(x) = e^{-x^2} \left[ \sin(x-1.4)+1 \right] \]

for $ x \in [-5,5] $ and then computes the value of $ \int x f(x)~dx / \int f(x)~dx $ using the Markov chain.

/* Example: ex_markov.cpp
-------------------------------------------------------------------
An example which demonstrates (i) the generation of an arbitrary
distribution through the creation of a Markov chain using the
Metropolis algorithm, and (ii) the evaluation of a integral over
that distribution using an expval object.
*/
#include <iostream>
#include <cmath>
#include <o2scl/test_mgr.h>
#include <o2scl/hist.h>
#include <o2scl/rng_gsl.h>
#include <o2scl/expval.h>
#include <o2scl/inte_qag_gsl.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 expsin_fun(double x) {
return exp(-x*x)*(sin(x-1.4)+1.0);
}
// Same function but now weighted with x
double expsin_fun2(double x) {
return x*exp(-x*x)*(sin(x-1.4)+1.0);
}
int main(int argc, char *argv[]) {
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
rng_gsl gr;
if (argc>=2) {
// We have to use the full o2scl::stoi() specification to
// avoid confusion with std::stoi()
gr.set_seed(o2scl::stoi(argv[1]));
}
// Construct histograms of the distribution. Initialize the grid
// spacing with a unform grid between -5 and 5.
hist h, h2;
// Compute the average value of the distribution, and use
// a expval_scalar object to estimate the uncertainty
static const size_t M=20;
expval_scalar sev(M,N*N/M);
// Store the coordinates to get the unblocked average value
std::vector<double> store;
// An initial point and the initial weight
double x_curr=0.01;
double w_curr=expsin_fun(x_curr);
double dev=0.0, dev2=0.0;
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_curr+r*step*2.0-step;
double w_new=expsin_fun(x_new);
// The Metrpolis algorithm
double r2=gr.random();
if (r2<w_new/w_curr && x_new>-5.0 && x_new<5.0) {
x_curr=x_new;
w_curr=w_new;
}
/*
There are two ways to do this, one is to just count up the
number of hits in each bin, the second is to record the
function value (the weight) in each bin, and divide by the
number of counts afterwards to get the average.
*/
// Update the histogram
h.update(x_curr);
// Update the second histogram with the weight
h2.update(x_curr,w_curr);
// Update the average value
sev.add(x_curr);
// Compute the average value using an simple unblocked average
store.push_back(x_curr);
}
}
// Normalize by the value of the distribution at x=1.05
double norm=h.get_wgt(1.05)/expsin_fun(1.05);
// Normalize the second histogram
for(size_t i=0;i<h.size();i++) {
if (h[i]>0.0) {
h2[i]/=h[i];
} else {
h2[i]=0.0;
}
}
// 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.setf(ios::left);
cout << i << " ";
cout.unsetf(ios::left);
cout.setf(ios::showpos);
cout << xx << " ";
cout.unsetf(ios::showpos);
cout << h[i]/norm << " ";
cout << h2[i] << " " << expsin_fun(xx) << endl;
}
dev+=fabs(h[i]/norm-expsin_fun(xx))/expsin_fun(xx);
dev2+=fabs(h2[i]-expsin_fun(xx))/expsin_fun(xx);
// Only test regions where we'll have enough statistics
if (expsin_fun(xx)>5.0e-3) {
t.test_rel(h.get_wgt_i(i)/norm,expsin_fun(xx),0.1,"dist");
}
}
cout << endl;
// Compare the deviations of the first and the second
// histogram
cout << "Avg. devs.: " << dev/((double)h.size()) << " "
<< dev2/((double)h.size()) << endl;
cout << endl;
cout.precision(4);
// Result from direct integration
funct11 ff=expsin_fun;
funct11 ff2=expsin_fun2;
double exact=qag.integ(ff2,-5.0,5.0)/qag.integ(ff,-5.0,5.0);
cout << "Exact result: " << exact << endl;
// Output the average value and its uncertainty
double avg, std_dev, avg_err;
sev.current_avg(avg,std_dev,avg_err);
cout << "Avg, err, |avg-exact| inc. autocorr.: "
<< avg << " " << avg_err << " " << fabs(avg-exact) << endl;
t.test_abs(avg,exact,3.0*avg_err,"Average value");
// Compute average and uncertainty ignoring autocorrelations
double avg2, std_dev2, avg_err2;
avg2=vector_mean(store.size(),store);
std_dev2=vector_variance(store.size(),store);
avg_err2=vector_variance(store.size(),store)/
sqrt(((double)store.size()));
cout << "Avg, err, |avg-exact| w/o autocorr. : "
<< avg2 << " " << avg_err2 << " " << fabs(avg2-exact) << endl;
// Show that expval and vector_mean() report the same average
t.test_rel(avg,avg2,1.0e-12,"Consistency");
cout << endl;
t.report();
return 0;
}
// End of example

Typical output is given below. Ignoring autocorrelations underestimates the uncertainty. The o2scl::expval_scalar class includes the effect of autocorrelations by forming block averages, and reports the correct uncertainty. (This is only the output from one random seed. This same ordering will hold for most but not all seeds.)

0  -4.950000e+00 0.000000e+00 0.000000e+00 2.131525e-11
3  -4.650000e+00 0.000000e+00 0.000000e+00 5.009023e-10
6  -4.350000e+00 0.000000e+00 0.000000e+00 9.131547e-09
9  -4.050000e+00 0.000000e+00 0.000000e+00 1.309343e-07
12 -3.750000e+00 4.095279e-06 1.879347e-06 1.488688e-06
15 -3.450000e+00 1.228584e-05 1.593305e-05 1.348287e-05
18 -3.150000e+00 5.323863e-05 9.847965e-05 9.747130e-05
21 -2.850000e+00 5.487674e-04 6.009592e-04 5.624059e-04
24 -2.550000e+00 2.706979e-03 2.637892e-03 2.584240e-03
27 -2.250000e+00 9.587048e-03 9.593383e-03 9.410936e-03
30 -1.950000e+00 2.736465e-02 2.725846e-02 2.693191e-02
33 -1.650000e+00 5.785401e-02 5.999921e-02 5.970011e-02
36 -1.350000e+00 1.004859e-01 1.000317e-01 9.993669e-02
39 -1.050000e+00 1.204053e-01 1.200373e-01 1.202766e-01
42 -7.500000e-01 9.158682e-02 9.303518e-02 9.293227e-02
45 -4.500000e-01 3.099307e-02 3.287070e-02 3.162602e-02
48 -1.500000e-01 6.757210e-04 1.473156e-03 2.114248e-04
51 +1.500000e-01 4.999517e-02 5.158937e-02 4.988035e-02
54 +4.500000e-01 1.519635e-01 1.529699e-01 1.523810e-01
57 +7.500000e-01 2.238029e-01 2.246116e-01 2.249580e-01
60 +1.050000e+00 2.181842e-01 2.179783e-01 2.181842e-01
63 +1.350000e+00 1.529259e-01 1.538486e-01 1.535435e-01
66 +1.650000e+00 8.190148e-02 8.253733e-02 8.196725e-02
69 +1.950000e+00 3.340519e-02 3.447211e-02 3.397864e-02
72 +2.250000e+00 1.076239e-02 1.132212e-02 1.108511e-02
75 +2.550000e+00 3.153365e-03 2.930076e-03 2.868544e-03
78 +2.850000e+00 7.576266e-04 6.258173e-04 5.914089e-04
81 +3.150000e+00 5.733391e-05 1.046556e-04 9.733109e-05
84 +3.450000e+00 4.095279e-06 1.731409e-05 1.278395e-05
87 +3.750000e+00 0.000000e+00 0.000000e+00 1.336916e-06
90 +4.050000e+00 0.000000e+00 0.000000e+00 1.107648e-07
93 +4.350000e+00 0.000000e+00 0.000000e+00 7.207155e-09
96 +4.650000e+00 0.000000e+00 0.000000e+00 3.628586e-10
99 +4.950000e+00 0.000000e+00 0.000000e+00 1.376924e-11

Avg. devs.: 4.197450e-01 3.578040e-01

Exact result: 2.8463e-01
Avg, err, |avg-exact| inc. autocorr.: 2.8739e-01 7.6799e-03 2.7640e-03
Avg, err, |avg-exact| w/o autocorr. : 2.8739e-01 1.2415e-03 2.7640e-03

48 tests performed.
All tests passed.

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..