pub fn integrate<T, F>(
f: &F,
a: &[f64],
b: &[f64],
gk_order: usize,
tci_options: TCI2Options,
) -> Result<T>Expand description
Integrate a function over a hypercube [a, b] using TCI and
Gauss-Kronrod quadrature.
The integrand is sampled at Gauss-Kronrod nodes in each dimension,
and a tensor-train approximation is built via crossinterpolate2.
The integral is then computed as a weighted sum of the tensor train.
§Arguments
f– function to integrate, takes a coordinate slice&[f64]a– lower bounds for each dimensionb– upper bounds for each dimensiongk_order– Gauss-Kronrod quadrature order (15 or 31)tci_options– options for the TCI approximation
§Returns
The approximate integral value.
§Errors
Returns an error if a and b have different lengths, or if
gk_order is not 15 or 31.
§Examples
use tensor4all_tensorci::integration::integrate;
use tensor4all_tensorci::TCI2Options;
// Integrate (x^2 + y^2) over [0,1]^2 = 2/3
let f = |x: &[f64]| x.iter().map(|&xi| xi * xi).sum::<f64>();
let a = vec![0.0, 0.0];
let b = vec![1.0, 1.0];
let options = TCI2Options { tolerance: 1e-10, ..TCI2Options::default() };
let result: f64 = integrate(&f, &a, &b, 15, options).unwrap();
assert!((result - 2.0 / 3.0).abs() < 1e-8);