use crate::Value;
use ffi::FFI;
use std::marker::PhantomData;
use std::mem::transmute;
use std::os::raw::c_void;
use std::slice;
ffi_wrapper!(PlainMonteCarlo, *mut sys::gsl_monte_plain_state, gsl_monte_plain_free,
"The plain Monte Carlo algorithm samples points randomly from the integration region to estimate
the integral and its error. Using this algorithm the estimate of the integral E(f; N) for N
randomly distributed points x_i is given by,
```text
E(f; N) = = V <f> = (V / N) sum_i^N f(x_i)
```
where V is the volume of the integration region. The error on this estimate `sigma(E;N)` is
calculated from the estimated variance of the mean,
```text
sigma^2 (E; N) = (V^2 / N^2) sum_i^N (f(x_i) - <f>)^2.
```
For large N this variance decreases asymptotically as `Var(f)/N`, where `Var(f)` is the true
variance of the function over the integration region. The error estimate itself should decrease as
`sigma(f)/sqrt{N}`. The familiar law of errors decreasing as `1/sqrt{N}` applies-to reduce the
error by a factor of 10 requires a 100-fold increase in the number of sample points.");
impl PlainMonteCarlo {
#[doc(alias = "gsl_monte_plain_alloc")]
pub fn new(dim: usize) -> Option<PlainMonteCarlo> {
let tmp = unsafe { sys::gsl_monte_plain_alloc(dim) };
if tmp.is_null() {
None
} else {
Some(PlainMonteCarlo::wrap(tmp))
}
}
#[doc(alias = "gsl_monte_plain_init")]
pub fn init(&mut self) -> Result<(), Value> {
let ret = unsafe { sys::gsl_monte_plain_init(self.unwrap_unique()) };
result_handler!(ret, ())
}
#[doc(alias = "gsl_monte_plain_integrate")]
pub fn integrate<F: FnMut(&[f64]) -> f64>(
&mut self,
f: F,
xl: &[f64],
xu: &[f64],
t_calls: usize,
r: &mut crate::Rng,
) -> Result<(f64, f64), Value> {
assert!(xl.len() == xu.len());
let mut result = 0f64;
let mut abserr = 0f64;
let f: Box<F> = Box::new(f);
let ret = unsafe {
let func = sys::gsl_monte_function {
f: transmute(monte_trampoline::<F> as usize),
dim: xl.len() as _,
params: Box::into_raw(f) as *mut _,
};
sys::gsl_monte_plain_integrate(
&func,
xl.as_ptr(),
xu.as_ptr(),
xl.len() as _,
t_calls,
r.unwrap_unique(),
self.unwrap_unique(),
&mut result,
&mut abserr,
)
};
result_handler!(ret, (result, abserr))
}
}
ffi_wrapper!(MiserMonteCarlo, *mut sys::gsl_monte_miser_state, gsl_monte_miser_free,
"The MISER algorithm of Press and Farrar is based on recursive stratified sampling. This technique
aims to reduce the overall integration error by concentrating integration points in the regions of
highest variance.
The idea of stratified sampling begins with the observation that for two disjoint regions a and b
with Monte Carlo estimates of the integral E_a(f) and E_b(f) and variances `sigma_a^2(f)` and
`sigma_b^2(f)`, the variance `Var(f)` of the combined estimate `E(f) = (1/2) (E_a(f) + E_b(f))`
is given by,
```text
Var(f) = (sigma_a^2(f) / 4 N_a) + (sigma_b^2(f) / 4 N_b).
```
It can be shown that this variance is minimized by distributing the points such that,
```text
N_a / (N_a + N_b) = sigma_a / (sigma_a + sigma_b).
```
Hence the smallest error estimate is obtained by allocating sample points in proportion to the
standard deviation of the function in each sub-region.
The MISER algorithm proceeds by bisecting the integration region along one coordinate axis to give
two sub-regions at each step. The direction is chosen by examining all d possible bisections and
selecting the one which will minimize the combined variance of the two sub-regions. The variance in
the sub-regions is estimated by sampling with a fraction of the total number of points available to
the current step. The same procedure is then repeated recursively for each of the two half-spaces
from the best bisection. The remaining sample points are allocated to the sub-regions using the
formula for N_a and N_b. This recursive allocation of integration points continues down to a
user-specified depth where each sub-region is integrated using a plain Monte Carlo estimate. These
individual values and their error estimates are then combined upwards to give an overall result and
an estimate of its error.");
impl MiserMonteCarlo {
#[doc(alias = "gsl_monte_miser_alloc")]
pub fn new(dim: usize) -> Option<MiserMonteCarlo> {
let tmp = unsafe { sys::gsl_monte_miser_alloc(dim) };
if tmp.is_null() {
None
} else {
Some(MiserMonteCarlo::wrap(tmp))
}
}
#[doc(alias = "gsl_monte_miser_init")]
pub fn init(&mut self) -> Result<(), Value> {
let ret = unsafe { sys::gsl_monte_miser_init(self.unwrap_unique()) };
result_handler!(ret, ())
}
#[doc(alias = "gsl_monte_miser_integrate")]
pub fn integrate<F: FnMut(&[f64]) -> f64>(
&mut self,
f: F,
xl: &[f64],
xu: &[f64],
t_calls: usize,
r: &mut crate::Rng,
) -> Result<(f64, f64), Value> {
assert!(xl.len() == xu.len());
let mut result = 0f64;
let mut abserr = 0f64;
let f: Box<F> = Box::new(f);
let ret = unsafe {
let mut func = sys::gsl_monte_function {
f: transmute(monte_trampoline::<F> as usize),
dim: xl.len() as _,
params: Box::into_raw(f) as *mut _,
};
sys::gsl_monte_miser_integrate(
&mut func,
xl.as_ptr(),
xu.as_ptr(),
xl.len() as _,
t_calls,
r.unwrap_unique(),
self.unwrap_unique(),
&mut result,
&mut abserr,
)
};
result_handler!(ret, (result, abserr))
}
#[doc(alias = "gsl_monte_miser_params_get")]
pub fn get_params(&self) -> MiserParams {
let mut m = sys::gsl_monte_miser_params {
estimate_frac: 0f64,
min_calls: 0,
min_calls_per_bisection: 0,
alpha: 0f64,
dither: 0f64,
};
unsafe {
sys::gsl_monte_miser_params_get(self.unwrap_shared(), &mut m);
}
MiserParams(m)
}
#[doc(alias = "gsl_monte_miser_params_set")]
pub fn set_params(&mut self, params: &MiserParams) {
unsafe {
sys::gsl_monte_miser_params_set(self.unwrap_unique(), ¶ms.0 as *const _);
}
}
}
#[derive(Debug, Clone)]
#[repr(C)]
pub struct MiserParams(pub sys::gsl_monte_miser_params);
ffi_wrapper!(VegasMonteCarlo, *mut sys::gsl_monte_vegas_state, gsl_monte_vegas_free,
"The VEGAS algorithm of Lepage is based on importance sampling. It samples points from the probability
distribution described by the function |f|, so that the points are concentrated in the regions that
make the largest contribution to the integral.
In general, if the Monte Carlo integral of f is sampled with points distributed according to a
probability distribution described by the function g, we obtain an estimate E_g(f; N),
```text
E_g(f; N) = E(f/g; N)
```
with a corresponding variance,
```text
Var_g(f; N) = Var(f/g; N).
```
If the probability distribution is chosen as g = |f|/I(|f|) then it can be shown that the variance
V_g(f; N) vanishes, and the error in the estimate will be zero. In practice it is not possible to
sample from the exact distribution g for an arbitrary function, so importance sampling algorithms
aim to produce efficient approximations to the desired distribution.
The VEGAS algorithm approximates the exact distribution by making a number of passes over the
integration region while histogramming the function f. Each histogram is used to define a sampling
distribution for the next pass. Asymptotically this procedure converges to the desired distribution.
In order to avoid the number of histogram bins growing like K^d the probability distribution is
approximated by a separable function: g(x_1, x_2, ...) = g_1(x_1) g_2(x_2) ... so that the number
of bins required is only Kd. This is equivalent to locating the peaks of the function from the
projections of the integrand onto the coordinate axes. The efficiency of VEGAS depends on the
validity of this assumption. It is most efficient when the peaks of the integrand are
well-localized. If an integrand can be rewritten in a form which is approximately separable this
will increase the efficiency of integration with VEGAS.
VEGAS incorporates a number of additional features, and combines both stratified sampling and
importance sampling. The integration region is divided into a number of “boxes”, with each box
getting a fixed number of points (the goal is 2). Each box can then have a fractional number of
bins, but if the ratio of bins-per-box is less than two, Vegas switches to a kind variance reduction
(rather than importance sampling).
The VEGAS algorithm computes a number of independent estimates of the integral internally, according
to the iterations parameter described below, and returns their weighted average. Random sampling of
the integrand can occasionally produce an estimate where the error is zero, particularly if the
function is constant in some regions. An estimate with zero error causes the weighted average to
break down and must be handled separately. In the original Fortran implementations of VEGAS the
error estimate is made non-zero by substituting a small value (typically 1e-30). The implementation
in GSL differs from this and avoids the use of an arbitrary constant—it either assigns the value a
weight which is the average weight of the preceding estimates or discards it according to the
following procedure,
current estimate has zero error, weighted average has finite error
* The current estimate is assigned a weight which is the average weight of the preceding estimates.
current estimate has finite error, previous estimates had zero error
* The previous estimates are discarded and the weighted averaging procedure begins with the current estimate.
current estimate has zero error, previous estimates had zero error
* The estimates are averaged using the arithmetic mean, but no error is computed.");
impl VegasMonteCarlo {
#[doc(alias = "gsl_monte_vegas_alloc")]
pub fn new(dim: usize) -> Option<VegasMonteCarlo> {
let tmp = unsafe { sys::gsl_monte_vegas_alloc(dim) };
if tmp.is_null() {
None
} else {
Some(VegasMonteCarlo::wrap(tmp))
}
}
#[doc(alias = "gsl_monte_vegas_init")]
pub fn init(&mut self) -> Result<(), Value> {
let ret = unsafe { sys::gsl_monte_vegas_init(self.unwrap_unique()) };
result_handler!(ret, ())
}
#[doc(alias = "gsl_monte_vegas_integrate")]
pub fn integrate<F: FnMut(&[f64]) -> f64>(
&mut self,
f: F,
xl: &[f64],
xu: &[f64],
t_calls: usize,
r: &mut crate::Rng,
) -> Result<(f64, f64), Value> {
assert!(xl.len() == xu.len());
let mut result = 0f64;
let mut abserr = 0f64;
let f: Box<F> = Box::new(f);
let ret = unsafe {
let mut func = sys::gsl_monte_function {
f: transmute(monte_trampoline::<F> as usize),
dim: xl.len() as _,
params: Box::into_raw(f) as *mut _,
};
sys::gsl_monte_vegas_integrate(
&mut func,
xl.as_ptr() as usize as *mut _,
xu.as_ptr() as usize as *mut _,
xl.len() as _,
t_calls,
r.unwrap_unique(),
self.unwrap_unique(),
&mut result,
&mut abserr,
)
};
result_handler!(ret, (result, abserr))
}
#[doc(alias = "gsl_monte_vegas_chisq")]
pub fn chisq(&mut self) -> f64 {
unsafe { sys::gsl_monte_vegas_chisq(self.unwrap_unique()) }
}
#[doc(alias = "gsl_monte_vegas_runval")]
pub fn runval(&mut self) -> (f64, f64) {
let mut result = 0.;
let mut sigma = 0.;
unsafe { sys::gsl_monte_vegas_runval(self.unwrap_unique(), &mut result, &mut sigma) };
(result, sigma)
}
#[doc(alias = "gsl_monte_vegas_params_get")]
pub fn get_params(&self) -> VegasParams {
let mut params = VegasParams::default();
unsafe {
sys::gsl_monte_vegas_params_get(self.unwrap_shared(), &mut params.inner as *mut _);
}
params
}
#[doc(alias = "gsl_monte_vegas_params_set")]
pub fn set_params(&mut self, params: &VegasParams) {
unsafe {
sys::gsl_monte_vegas_params_set(self.unwrap_unique(), ¶ms.inner as *const _);
}
}
}
pub struct VegasParams<'a> {
inner: sys::gsl_monte_vegas_params,
lt: PhantomData<&'a ()>,
}
impl<'a> VegasParams<'a> {
pub fn new(
alpha: f64,
iterations: usize,
stage: i32,
mode: ::VegasMode,
verbosity: VegasVerbosity,
stream: Option<&'a mut ::IOStream>,
) -> Result<VegasParams, String> {
if !verbosity.is_off() && stream.is_none() {
return Err(
"rust-GSL: need to provide an input stream for Vegas Monte Carlo \
integration if verbosity is not 'Off'"
.to_string(),
);
} else if verbosity.is_off() && stream.is_some() {
return Err(
"rust-GSL: need to provide the verbosity flag for Vegas Monta Carlo \
integration, currently set to 'Off'"
.to_string(),
);
}
let stream = if let Some(stream) = stream {
if !stream.write_mode() {
return Err("rust-GSL: input stream not flagged as 'write' mode".to_string());
}
stream.as_raw()
} else {
std::ptr::null_mut()
};
Ok(VegasParams {
inner: sys::gsl_monte_vegas_params {
alpha,
iterations,
stage,
mode: mode.into(),
verbose: verbosity.to_int(),
ostream: stream,
},
lt: PhantomData,
})
}
}
impl<'a> std::default::Default for VegasParams<'a> {
fn default() -> VegasParams<'a> {
VegasParams {
inner: sys::gsl_monte_vegas_params {
alpha: 1.5,
iterations: 5,
stage: 0,
mode: ::VegasMode::ImportanceOnly.into(),
verbose: -1,
ostream: std::ptr::null_mut(),
},
lt: PhantomData,
}
}
}
#[derive(Clone, Copy)]
pub enum VegasVerbosity {
Off, Summary, Grid, Rebinning, }
impl VegasVerbosity {
fn to_int(self) -> i32 {
match self {
VegasVerbosity::Off => -1,
VegasVerbosity::Summary => 0,
VegasVerbosity::Grid => 1,
VegasVerbosity::Rebinning => 2,
}
}
fn is_off(self) -> bool {
matches!(self, VegasVerbosity::Off)
}
}
unsafe extern "C" fn monte_trampoline<F: FnMut(&[f64]) -> f64>(
x: *mut f64,
dim: usize,
param: *mut c_void,
) -> f64 {
let f: &mut F = &mut *(param as *mut F);
f(slice::from_raw_parts(x, dim))
}
#[test]
fn plain() {
use std::f64::consts::PI;
fn g(k: &[f64]) -> f64 {
let a = 1f64 / (PI * PI * PI);
a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
}
let xl: [f64; 3] = [0f64; 3];
let xu: [f64; 3] = [PI, PI, PI];
let calls = 500000;
crate::RngType::env_setup();
let mut r = crate::Rng::new(::RngType::default()).unwrap();
{
let mut s = PlainMonteCarlo::new(3).unwrap();
let (res, err) = s.integrate(g, &xl, &xu, calls, &mut r).unwrap();
assert_eq!(&format!("{:.6}", res), "1.412209");
assert_eq!(&format!("{:.6}", err), "0.013436");
}
}
#[test]
fn miser() {
use std::f64::consts::PI;
fn g(k: &[f64]) -> f64 {
let a = 1f64 / (PI * PI * PI);
a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
}
let xl: [f64; 3] = [0f64; 3];
let xu: [f64; 3] = [PI, PI, PI];
let calls = 500000;
crate::RngType::env_setup();
let mut r = crate::Rng::new(::RngType::default()).unwrap();
{
let mut s = MiserMonteCarlo::new(3).unwrap();
let (res, err) = s.integrate(g, &xl, &xu, calls, &mut r).unwrap();
assert_eq!(&format!("{:.6}", res), "1.389530");
assert_eq!(&format!("{:.6}", err), "0.005011");
}
}
#[test]
fn miser_closure() {
use std::f64::consts::PI;
let xl: [f64; 3] = [0f64; 3];
let xu: [f64; 3] = [PI, PI, PI];
let calls = 500000;
crate::RngType::env_setup();
let mut r = crate::Rng::new(::RngType::default()).unwrap();
{
let mut s = MiserMonteCarlo::new(3).unwrap();
let (res, err) = s
.integrate(
|k| {
let a = 1f64 / (PI * PI * PI);
a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
},
&xl,
&xu,
calls,
&mut r,
)
.unwrap();
assert_eq!(&format!("{:.6}", res), "1.389530");
assert_eq!(&format!("{:.6}", err), "0.005011");
}
}
#[test]
fn vegas_warm_up() {
use std::f64::consts::PI;
fn g(k: &[f64]) -> f64 {
let a = 1f64 / (PI * PI * PI);
a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
}
let xl: [f64; 3] = [0f64; 3];
let xu: [f64; 3] = [PI, PI, PI];
crate::RngType::env_setup();
let mut r = crate::Rng::new(::RngType::default()).unwrap();
{
let mut s = VegasMonteCarlo::new(3).unwrap();
let (res, err) = s.integrate(g, &xl, &xu, 10000, &mut r).unwrap();
assert_eq!(&format!("{:.6}", res), "1.385603");
assert_eq!(&format!("{:.6}", err), "0.002212");
}
}
#[test]
fn vegas() {
use std::f64::consts::PI;
fn g(k: &[f64]) -> f64 {
let a = 1. / (PI * PI * PI);
a / (1. - k[0].cos() * k[1].cos() * k[2].cos())
}
let calls = 500000;
let xl: [f64; 3] = [0f64; 3];
let xu: [f64; 3] = [PI, PI, PI];
crate::RngType::env_setup();
let mut r = crate::Rng::new(::RngType::default()).unwrap();
{
let mut s = VegasMonteCarlo::new(3).unwrap();
s.integrate(g, &xl, &xu, 10000, &mut r).unwrap();
let mut res;
let mut err;
loop {
let (_res, _err) = s.integrate(g, &xl, &xu, calls / 5, &mut r).unwrap();
res = _res;
err = _err;
println!(
"result = {:.6} sigma = {:.6} chisq/dof = {:.1}",
res,
err,
s.chisq()
);
if (s.chisq() - 1f64).abs() <= 0.5f64 {
break;
}
}
assert_eq!(&format!("{:.6}", res), "1.393307");
assert_eq!(&format!("{:.6}", err), "0.000335");
}
}