Expand description
§bilby
A high-performance numerical quadrature (integration) library for Rust.
bilby provides Gaussian quadrature rules, adaptive integration, and multi-dimensional cubature methods with a consistent, generic API.
§Design
- Generic over
F: Floatvianum_traits::Float— works withf32,f64, and compatible types (e.g., echidna AD types) - Precomputed rules — generate nodes and weights once, integrate many times
- No heap allocation on hot paths — rule construction allocates, integration does not
- Separate 1D and N-D APIs —
Fn(F) -> Ffor 1D,Fn(&[F]) -> Ffor N-D no_stdcompatible — works without the standard library (withalloc)
§Feature Flags
| Feature | Default | Description |
|---|---|---|
std | Yes | Enables std::error::Error impl and cache module |
parallel | No | Enables rayon-based _par methods (implies std) |
§Quick Start
use bilby::GaussLegendre;
// Create a 10-point Gauss-Legendre rule
let gl = GaussLegendre::new(10).unwrap();
// Integrate x^2 over [0, 1] (exact result = 1/3)
let result = gl.rule().integrate(0.0, 1.0, |x: f64| x * x);
assert!((result - 1.0 / 3.0).abs() < 1e-14);§Gauss-Kronrod Error Estimation
use bilby::{GaussKronrod, GKPair};
let gk = GaussKronrod::new(GKPair::G7K15);
let (estimate, error) = gk.integrate(0.0, core::f64::consts::PI, f64::sin);
assert!((estimate - 2.0).abs() < 1e-14);§Adaptive Integration
use bilby::adaptive_integrate;
let result = adaptive_integrate(|x: f64| x.sin(), 0.0, core::f64::consts::PI, 1e-12).unwrap();
assert!((result.value - 2.0).abs() < 1e-12);
assert!(result.is_converged());§Infinite Domains
use bilby::integrate_infinite;
// Integral of exp(-x^2) over (-inf, inf) = sqrt(pi)
let result = integrate_infinite(|x: f64| (-x * x).exp(), 1e-10).unwrap();
assert!((result.value - core::f64::consts::PI.sqrt()).abs() < 1e-8);§Multi-Dimensional Integration
use bilby::cubature::{TensorProductRule, adaptive_cubature};
use bilby::GaussLegendre;
// Tensor product: 10-point GL in each of 2 dimensions
let gl = GaussLegendre::new(10).unwrap();
let tp = TensorProductRule::isotropic(gl.rule(), 2).unwrap();
let result = tp.rule().integrate_box(
&[0.0, 0.0], &[1.0, 1.0],
|x| x[0] * x[1],
);
assert!((result - 0.25).abs() < 1e-14);§Precomputed Rule Cache
Available when the std feature is enabled (default):
use bilby::cache::GL10;
let result = GL10.rule().integrate(0.0, 1.0, |x: f64| x * x);
assert!((result - 1.0 / 3.0).abs() < 1e-14);§Tanh-Sinh (Double Exponential)
use bilby::tanh_sinh_integrate;
// Endpoint singularity: integral of 1/sqrt(x) over [0, 1] = 2
let result = tanh_sinh_integrate(|x| 1.0 / x.sqrt(), 0.0, 1.0, 1e-10).unwrap();
assert!((result.value - 2.0).abs() < 1e-7);§Cauchy Principal Value
use bilby::pv_integrate;
// PV integral of x^2/(x - 0.3) over [0, 1]
let exact = 0.8 + 0.09 * (7.0_f64 / 3.0).ln();
let result = pv_integrate(|x| x * x, 0.0, 1.0, 0.3, 1e-10).unwrap();
assert!((result.value - exact).abs() < 1e-7);§Oscillatory Integration
use bilby::integrate_oscillatory_sin;
// Integral of sin(100x) over [0, 1]
let exact = (1.0 - 100.0_f64.cos()) / 100.0;
let result = integrate_oscillatory_sin(|_| 1.0, 0.0, 1.0, 100.0, 1e-10).unwrap();
assert!((result.value - exact).abs() < 1e-8);§Weighted Integration
use bilby::weighted::{weighted_integrate, WeightFunction};
// Gauss-Hermite: integral of e^(-x^2) over (-inf, inf) = sqrt(pi)
let result = weighted_integrate(|_| 1.0, WeightFunction::Hermite, 20).unwrap();
assert!((result - core::f64::consts::PI.sqrt()).abs() < 1e-12);§Sparse Grids
use bilby::cubature::SparseGrid;
let sg = SparseGrid::clenshaw_curtis(3, 3).unwrap();
let result = sg.rule().integrate_box(
&[0.0, 0.0, 0.0], &[1.0, 1.0, 1.0],
|x| x[0] * x[1] * x[2],
);
assert!((result - 0.125).abs() < 1e-12);Re-exports§
pub use adaptive::adaptive_integrate;pub use adaptive::adaptive_integrate_with_breaks;pub use adaptive::AdaptiveIntegrator;pub use cauchy_pv::pv_integrate;pub use cauchy_pv::CauchyPV;pub use clenshaw_curtis::ClenshawCurtis;pub use error::QuadratureError;pub use gauss_chebyshev::GaussChebyshevFirstKind;pub use gauss_chebyshev::GaussChebyshevSecondKind;pub use gauss_hermite::GaussHermite;pub use gauss_jacobi::GaussJacobi;pub use gauss_kronrod::GKPair;pub use gauss_kronrod::GaussKronrod;pub use gauss_laguerre::GaussLaguerre;pub use gauss_legendre::GaussLegendre;pub use gauss_lobatto::GaussLobatto;pub use gauss_radau::GaussRadau;pub use oscillatory::integrate_oscillatory_cos;pub use oscillatory::integrate_oscillatory_sin;pub use oscillatory::OscillatoryIntegrator;pub use oscillatory::OscillatoryKernel;pub use result::QuadratureResult;pub use rule::QuadratureRule;pub use tanh_sinh::tanh_sinh_integrate;pub use tanh_sinh::TanhSinh;pub use transforms::integrate_infinite;pub use transforms::integrate_semi_infinite_lower;pub use transforms::integrate_semi_infinite_upper;pub use weighted::weighted_integrate;pub use weighted::WeightFunction;pub use weighted::WeightedIntegrator;
Modules§
- adaptive
- Adaptive quadrature via global subdivision.
- cache
- Precomputed quadrature rule cache.
- cauchy_
pv - Cauchy principal value integration.
- clenshaw_
curtis - Clenshaw-Curtis quadrature rule.
- cubature
- Multi-dimensional integration (cubature).
- error
- Error types for quadrature operations.
- gauss_
chebyshev - Gauss-Chebyshev quadrature rules (Types I and II).
- gauss_
hermite - Gauss-Hermite quadrature rule.
- gauss_
jacobi - Gauss-Jacobi quadrature rule.
- gauss_
kronrod - Gauss-Kronrod embedded quadrature pairs.
- gauss_
laguerre - Gauss-Laguerre quadrature rule.
- gauss_
legendre - Gauss-Legendre quadrature rule.
- gauss_
lobatto - Gauss-Lobatto quadrature rule.
- gauss_
radau - Gauss-Radau quadrature rule.
- oscillatory
- Oscillatory integration via Filon-Clenshaw-Curtis.
- result
- Result types for quadrature operations.
- rule
- Core quadrature rule types.
- tanh_
sinh - Tanh-sinh (double-exponential) quadrature.
- transforms
- Domain transforms for semi-infinite and infinite intervals.
- weighted
- Weighted integration: ∫ f(x) · w(x) dx for known weight functions.