gauss_quad/
lib.rs

1//! # gauss-quad
2//!
3//! **gauss-quad** is a [Gaussian quadrature](https://en.wikipedia.org/wiki/Gaussian_quadrature) library for numerical integration.
4//!
5//! ## Quadrature rules
6//!
7//! **gauss-quad** implements the following quadrature rules:
8//! * [Gauss-Legendre](https://en.wikipedia.org/wiki/Gauss%E2%80%93Legendre_quadrature)
9//! * [Gauss-Jacobi](https://en.wikipedia.org/wiki/Gauss%E2%80%93Jacobi_quadrature)
10//! * [Gauss-Laguerre](https://en.wikipedia.org/wiki/Gauss%E2%80%93Laguerre_quadrature) (generalized)
11//! * [Gauss-Hermite](https://en.wikipedia.org/wiki/Gauss%E2%80%93Hermite_quadrature)
12//! * [Gauss-Chebyshev](https://en.wikipedia.org/wiki/Chebyshev%E2%80%93Gauss_quadrature)
13//! * [Midpoint](https://en.wikipedia.org/wiki/Riemann_sum#Midpoint_rule)
14//! * [Simpson](https://en.wikipedia.org/wiki/Simpson%27s_rule)
15//!
16//! ## Using **gauss-quad**
17//!
18//! To use any of the quadrature rules in your project, first initialize the rule with
19//! a specified degree and then you can use it for integration, e.g.:
20//! ```
21//! use gauss_quad::GaussLegendre;
22//! # use gauss_quad::legendre::GaussLegendreError;
23//! // This macro is used in these docs to compare floats.
24//! // The assertion succeeds if the two sides are within floating point error,
25//! // or an optional epsilon.
26//! use approx::assert_abs_diff_eq;
27//!
28//! // initialize the quadrature rule
29//! let degree = 10;
30//! let quad = GaussLegendre::new(degree)?;
31//!
32//! // use the rule to integrate a function
33//! let left_bound = 0.0;
34//! let right_bound = 1.0;
35//! let integral = quad.integrate(left_bound, right_bound, |x| x * x);
36//! assert_abs_diff_eq!(integral, 1.0 / 3.0);
37//! # Ok::<(), GaussLegendreError>(())
38//! ```
39//!
40//! ## Setting up a quadrature rule
41//!
42//! Using a quadrature rule takes two steps:
43//! 1. Initialization
44//! 2. Integration
45//!
46//! First, rules must be initialized using some specific input parameters.
47//!
48//! Then, you can integrate functions using those rules:
49//! ```
50//! # use gauss_quad::*;
51//! # let degree = 5;
52//! # let alpha = 1.2;
53//! # let beta = 1.2;
54//! # let a = 0.0;
55//! # let b = 1.0;
56//! # let c = -10.;
57//! # let d = 100.;
58//! let gauss_legendre = GaussLegendre::new(degree)?;
59//! // Integrate on the domain [a, b]
60//! let x_cubed = gauss_legendre.integrate(a, b, |x| x * x * x);
61//!
62//! let gauss_jacobi = GaussJacobi::new(degree, alpha, beta)?;
63//! // Integrate on the domain [c, d]
64//! let double_x = gauss_jacobi.integrate(c, d, |x| 2.0 * x);
65//!
66//! let gauss_laguerre = GaussLaguerre::new(degree, alpha)?;
67//! // no explicit domain, Gauss-Laguerre integration is done on the domain [0, ∞).
68//! let piecewise = gauss_laguerre.integrate(|x| if x > 0.0 && x < 2.0 { x } else { 0.0 });
69//!
70//! let gauss_hermite = GaussHermite::new(degree)?;
71//! // again, no explicit domain since Gauss-Hermite integration is done over the domain (-∞, ∞).
72//! let golden_polynomial = gauss_hermite.integrate(|x| x * x - x - 1.0);
73//! # Ok::<(), Box<dyn std::error::Error>>(())
74//! ```
75//!
76//! ## Specific quadrature rules
77//!
78//! Different rules may take different parameters.
79//!
80//! For example, the `GaussLaguerre` rule requires both a `degree` and an `alpha`
81//! parameter.
82//!
83//! `GaussLaguerre` is also defined as an improper integral over the domain [0, ∞).
84//! This means no domain bounds are needed in the `integrate` call.
85//! ```
86//! # use gauss_quad::laguerre::{GaussLaguerre, GaussLaguerreError};
87//! # use approx::assert_abs_diff_eq;
88//! # use core::f64::consts::PI;
89//! // initialize the quadrature rule
90//! let degree = 2;
91//! let alpha = 0.5;
92//! let quad = GaussLaguerre::new(degree, alpha)?;
93//!
94//! // use the rule to integrate a function
95//! let integral = quad.integrate(|x| x * x);
96//!
97//! assert_abs_diff_eq!(integral, 15.0 * PI.sqrt() / 8.0, epsilon = 1e-14);
98//! # Ok::<(), GaussLaguerreError>(())
99//! ```
100//!
101//! ## Errors
102//!
103//! Quadrature rules are only defined for a certain set of input values.
104//! For example, every rule is only defined for degrees where `degree > 1`.
105//! ```
106//! # use gauss_quad::GaussLaguerre;
107//! let degree = 1;
108//! assert!(GaussLaguerre::new(degree, 0.1).is_err());
109//! ```
110//!
111//! Specific rules may have other requirements.
112//! `GaussJacobi` for example, requires alpha and beta parameters larger than -1.0.
113//! ```
114//! # use gauss_quad::jacobi::GaussJacobi;
115//! let degree = 10;
116//! let alpha = 0.1;
117//! let beta = -1.1;
118//!
119//! assert!(GaussJacobi::new(degree, alpha, beta).is_err())
120//! ```
121//! Make sure to read the specific quadrature rule's documentation before using it.
122//!
123//! ## Passing functions to quadrature rules
124//!
125//! The `integrate` method takes any integrand that implements the [`Fn(f64) -> f64`](Fn) trait, i.e. functions of
126//! one `f64` parameter.
127//!
128//! ```
129//! # use gauss_quad::legendre::{GaussLegendre, GaussLegendreError};
130//! # use approx::assert_abs_diff_eq;
131//!
132//! // initialize the quadrature rule
133//! let degree = 2;
134//! let quad = GaussLegendre::new(degree)?;
135//!
136//! // use the rule to integrate a function
137//! let left_bound = 0.0;
138//! let right_bound = 1.0;
139//!
140//! let integral = quad.integrate(left_bound, right_bound, |x| x * x);
141//!
142//! assert_abs_diff_eq!(integral, 1.0 / 3.0);
143//! # Ok::<(), GaussLegendreError>(())
144//! ```
145//!
146//! ## Double integrals
147//!
148//! It is possible to use this crate to do double and higher integrals:
149//! ```
150//! # use gauss_quad::legendre::{GaussLegendre, GaussLegendreError};
151//! # use approx::assert_abs_diff_eq;
152//! let rule = GaussLegendre::new(3)?;
153//!
154//! // integrate x^2*y over the triangle in the xy-plane where x ϵ [0, 1] and y ϵ [0, x]:
155//! let double_int = rule.integrate(0.0, 1.0, |x| rule.integrate(0.0, x, |y| x * x * y));
156//!
157//! assert_abs_diff_eq!(double_int, 0.1);
158//! # Ok::<(), GaussLegendreError>(())
159//! ```
160//! However, the time complexity of the integration then scales with the number of nodes to
161//! the power of the depth of the integral, e.g. O(n³) for triple integrals.
162//!
163//! ## Feature flags
164//!
165//! `serde`: implements the [`Serialize`](serde::Serialize) and [`Deserialize`](serde::Deserialize) traits from
166//! the [`serde`](https://crates.io/crates/serde) crate for the quadrature rule structs.
167//!
168//! `rayon`: enables a parallel version of the `integrate` function on the quadrature rule structs. Can speed up integration if evaluating the integrand is expensive (takes ≫100 µs).
169
170// Only enable the nighlty `doc_auto_cfg` feature when
171// the `docsrs` configuration attribute is defined.
172// The config in Cargo.toml means that this is defined when we are on docs.rs (which uses the nightly compiler)
173// or if the environment variable RUSTDOCFLAGS is set as `RUSTDOCFLAGS="--cfg docsrs"`.
174// This lets us get a tag on docs.rs that says "Available on crate feature ... only".
175#![cfg_attr(docsrs, feature(doc_auto_cfg))]
176
177use nalgebra::{Dyn, Matrix, VecStorage};
178
179type DMatrixf64 = Matrix<f64, Dyn, Dyn, VecStorage<f64, Dyn, Dyn>>;
180
181pub mod chebyshev;
182mod data_api;
183mod gamma;
184pub mod hermite;
185pub mod jacobi;
186pub mod laguerre;
187pub mod legendre;
188pub mod midpoint;
189pub mod simpson;
190
191#[doc(inline)]
192pub use chebyshev::{GaussChebyshevFirstKind, GaussChebyshevSecondKind};
193#[doc(inline)]
194pub use data_api::{Node, Weight};
195#[doc(inline)]
196pub use hermite::GaussHermite;
197#[doc(inline)]
198pub use jacobi::GaussJacobi;
199#[doc(inline)]
200pub use laguerre::GaussLaguerre;
201#[doc(inline)]
202pub use legendre::GaussLegendre;
203#[doc(inline)]
204pub use midpoint::Midpoint;
205#[doc(inline)]
206pub use simpson::Simpson;