Struct gauss_quad::GaussLaguerre
source · pub struct GaussLaguerre {
pub nodes: Vec<f64>,
pub weights: Vec<f64>,
}
Expand description
A Gauss-Laguerre quadrature scheme.
These rules can perform integrals with integrands of the form x^alpha * e^(-x) * f(x) over the domain [0, ∞).
Example
Compute the factorial of 5:
// initialize a Gauss-Laguerre rule with 10 nodes
let quad = GaussLaguerre::init(10, 0.0);
// numerically evaluate this integral,
// which is a definition of the gamma function
let fact_5 = quad.integrate(|x| x.powi(5));
assert_abs_diff_eq!(fact_5, 1.0 * 2.0 * 3.0 * 4.0 * 5.0, epsilon = 1e-11);
Fields§
§nodes: Vec<f64>
§weights: Vec<f64>
Implementations§
source§impl GaussLaguerre
impl GaussLaguerre
sourcepub fn init(deg: usize, alpha: f64) -> GaussLaguerre
pub fn init(deg: usize, alpha: f64) -> GaussLaguerre
Initializes Gauss-Laguerre quadrature rule of the given degree by computing the nodes and weights
needed for the given alpha
parameter.
Panics
Panics if degree of quadrature is smaller than 2, or if alpha is smaller than -1
sourcepub fn nodes_and_weights(deg: usize, alpha: f64) -> (Vec<f64>, Vec<f64>)
pub fn nodes_and_weights(deg: usize, alpha: f64) -> (Vec<f64>, Vec<f64>)
Apply Golub-Welsch algorithm to determine Gauss-Laguerre nodes & weights construct companion matrix A for the Laguerre Polynomial using the relation: -n L_{n-1} + (2n+1) L_{n} -(n+1) L_{n+1} = x L_n The constructed matrix is symmetric and tridiagonal with (2n+1) on the diagonal & -(n+1) on the off-diagonal (n = row number). Root & weight finding are equivalent to eigenvalue problem. see Gil, Segura, Temme - Numerical Methods for Special Functions
Panics
Panics if degree of quadrature is smaller than 2, or if alpha is smaller than -1
Trait Implementations§
source§impl Clone for GaussLaguerre
impl Clone for GaussLaguerre
source§fn clone(&self) -> GaussLaguerre
fn clone(&self) -> GaussLaguerre
1.0.0 · source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source
. Read moresource§impl Debug for GaussLaguerre
impl Debug for GaussLaguerre
source§impl PartialEq for GaussLaguerre
impl PartialEq for GaussLaguerre
source§fn eq(&self, other: &GaussLaguerre) -> bool
fn eq(&self, other: &GaussLaguerre) -> bool
self
and other
values to be equal, and is used
by ==
.impl StructuralPartialEq for GaussLaguerre
Auto Trait Implementations§
impl RefUnwindSafe for GaussLaguerre
impl Send for GaussLaguerre
impl Sync for GaussLaguerre
impl Unpin for GaussLaguerre
impl UnwindSafe for GaussLaguerre
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
§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere SS: SubsetOf<SP>,
§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read more§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).§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.