Struct gauss_quad::GaussJacobi
source · pub struct GaussJacobi {
pub nodes: Vec<f64>,
pub weights: Vec<f64>,
}
Expand description
A Gauss-Jacobi quadrature scheme.
These rules can approximate integrals with singularities at the end points of the domain, [a, b].
Examples
// initialize the quadrature rule.
let quad = GaussJacobi::init(10, -0.5, 0.0);
// numerically integrate e^-x / sqrt(1 + x).
let integral = quad.integrate(-1.0, 1.0, |x| (-x).exp());
let dawson_function_of_sqrt_2 = 0.4525399074037225;
assert_abs_diff_eq!(integral, 2.0 * E * dawson_function_of_sqrt_2, epsilon = 1e-14);
Fields§
§nodes: Vec<f64>
§weights: Vec<f64>
Implementations§
source§impl GaussJacobi
impl GaussJacobi
sourcepub fn init(deg: usize, alpha: f64, beta: f64) -> GaussJacobi
pub fn init(deg: usize, alpha: f64, beta: f64) -> GaussJacobi
Initializes Gauss-Jacobi quadrature rule of the given degree by computing the nodes and weights
needed for the given alpha
and beta
.
Panics
Panics if degree of quadrature is smaller than 2, or if alpha or beta are smaller than -1
sourcepub fn nodes_and_weights(
deg: usize,
alpha: f64,
beta: f64
) -> (Vec<f64>, Vec<f64>)
pub fn nodes_and_weights( deg: usize, alpha: f64, beta: f64 ) -> (Vec<f64>, Vec<f64>)
Apply Golub-Welsch algorithm to determine Gauss-Jacobi nodes & weights see Gil, Segura, Temme - Numerical Methods for Special Functions
Panics
Panics if degree of quadrature is smaller than 2, or if alpha or beta are smaller than -1
Trait Implementations§
source§impl Clone for GaussJacobi
impl Clone for GaussJacobi
source§fn clone(&self) -> GaussJacobi
fn clone(&self) -> GaussJacobi
Returns a copy of the value. Read more
1.0.0 · source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
Performs copy-assignment from
source
. Read moresource§impl Debug for GaussJacobi
impl Debug for GaussJacobi
source§impl PartialEq for GaussJacobi
impl PartialEq for GaussJacobi
source§fn eq(&self, other: &GaussJacobi) -> bool
fn eq(&self, other: &GaussJacobi) -> bool
This method tests for
self
and other
values to be equal, and is used
by ==
.impl StructuralPartialEq for GaussJacobi
Auto Trait Implementations§
impl RefUnwindSafe for GaussJacobi
impl Send for GaussJacobi
impl Sync for GaussJacobi
impl Unpin for GaussJacobi
impl UnwindSafe for GaussJacobi
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
Mutably borrows from an owned value. Read more
§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>
The inverse inclusion map: attempts to construct
self
from the equivalent element of its
superset. Read more§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
Checks if
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
Use with care! Same as
self.to_subset
but without any property checks. Always succeeds.§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self
to the equivalent element of its superset.