pub struct GaussLegendre { /* private fields */ }
Expand description
A Gauss-Legendre quadrature scheme.
These rules can integrate functions on the domain [a, b].
§Examples
Basic usage:
// initialize a Gauss-Legendre rule with 2 nodes
let quad = GaussLegendre::new(2)?;
// numerically integrate x^2 - 1/3 over the domain [0, 1]
let integral = quad.integrate(0.0, 1.0, |x| x * x - 1.0 / 3.0);
assert_abs_diff_eq!(integral, 0.0);
The nodes and weights are computed in O(n) time, so large quadrature rules are feasible:
let quad = GaussLegendre::new(1_000_000)?;
let integral = quad.integrate(-3.0, 3.0, |x| x.sin());
assert_abs_diff_eq!(integral, 0.0);
Implementations§
Source§impl GaussLegendre
impl GaussLegendre
Sourcepub fn new(deg: usize) -> Result<Self, GaussLegendreError>
pub fn new(deg: usize) -> Result<Self, GaussLegendreError>
Initializes a Gauss-Legendre quadrature rule of the given degree by computing the needed nodes and weights.
A rule of degree n can integrate polynomials of degree 2n-1 exactly.
Uses the algorithm by Ignace Bogaert, which has linear time complexity.
§Errors
Returns an error if deg
is smaller than 2.
Sourcepub fn par_new(deg: usize) -> Result<Self, GaussLegendreError>
Available on crate feature rayon
only.
pub fn par_new(deg: usize) -> Result<Self, GaussLegendreError>
rayon
only.Sourcepub fn integrate<F>(&self, a: f64, b: f64, integrand: F) -> f64
pub fn integrate<F>(&self, a: f64, b: f64, integrand: F) -> f64
Perform quadrature integration of given integrand from a
to b
.
§Example
Basic usage
let glq_rule = GaussLegendre::new(3)?;
assert_abs_diff_eq!(glq_rule.integrate(-1.0, 1.0, |x| x.powi(5)), 0.0);
Sourcepub fn par_integrate<F>(&self, a: f64, b: f64, integrand: F) -> f64
Available on crate feature rayon
only.
pub fn par_integrate<F>(&self, a: f64, b: f64, integrand: F) -> f64
rayon
only.Source§impl GaussLegendre
impl GaussLegendre
Sourcepub fn nodes(&self) -> GaussLegendreNodes<'_> ⓘ
pub fn nodes(&self) -> GaussLegendreNodes<'_> ⓘ
Returns an iterator over the nodes of the quadrature rule.
Sourcepub fn weights(&self) -> GaussLegendreWeights<'_> ⓘ
pub fn weights(&self) -> GaussLegendreWeights<'_> ⓘ
Returns an iterator over the weights of the quadrature rule.
Sourcepub fn iter(&self) -> GaussLegendreIter<'_> ⓘ
pub fn iter(&self) -> GaussLegendreIter<'_> ⓘ
Returns an iterator over the node-weight pairs of the quadrature rule.
Sourcepub fn as_node_weight_pairs(&self) -> &[(Node, Weight)]
pub fn as_node_weight_pairs(&self) -> &[(Node, Weight)]
Returns a slice of all the node-weight pairs of the quadrature rule.
Sourcepub fn into_node_weight_pairs(self) -> Vec<(Node, Weight)>
pub fn into_node_weight_pairs(self) -> Vec<(Node, Weight)>
Converts the quadrature rule into a vector of node-weight pairs.
This function just returns the underlying vector without any computation or cloning.
Trait Implementations§
Source§impl Clone for GaussLegendre
impl Clone for GaussLegendre
Source§fn clone(&self) -> GaussLegendre
fn clone(&self) -> GaussLegendre
1.0.0 · Source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source
. Read moreSource§impl Debug for GaussLegendre
impl Debug for GaussLegendre
Source§impl<'de> Deserialize<'de> for GaussLegendre
impl<'de> Deserialize<'de> for GaussLegendre
Source§fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
Source§impl From<GaussLegendre> for GaussJacobi
Gauss-Legendre quadrature is equivalent to Gauss-Jacobi quadrature with alpha
= beta
= 0.
impl From<GaussLegendre> for GaussJacobi
Gauss-Legendre quadrature is equivalent to Gauss-Jacobi quadrature with alpha
= beta
= 0.
Source§fn from(value: GaussLegendre) -> Self
fn from(value: GaussLegendre) -> Self
Source§impl<'a> IntoIterator for &'a GaussLegendre
impl<'a> IntoIterator for &'a GaussLegendre
Source§impl IntoIterator for GaussLegendre
impl IntoIterator for GaussLegendre
Source§impl PartialEq for GaussLegendre
impl PartialEq for GaussLegendre
Source§impl Serialize for GaussLegendre
impl Serialize for GaussLegendre
impl StructuralPartialEq for GaussLegendre
Auto Trait Implementations§
impl Freeze for GaussLegendre
impl RefUnwindSafe for GaussLegendre
impl Send for GaussLegendre
impl Sync for GaussLegendre
impl Unpin for GaussLegendre
impl UnwindSafe for GaussLegendre
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left
is true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left(&self)
returns true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§impl<T> Pointable for T
impl<T> Pointable for T
Source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
Source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self
is actually part of its subset T
(and can be converted to it).Source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.