Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Gauss-Legendre integration.
Low-level mathematical utility that applies the Gauss-Legendre method to compute numeric integrals.
Derived types and interfaces | |
interface | sll_o_gauss_legendre_integrate_1d |
Compute numerical integration with Gauss-Legendre formula. More... | |
Functions/Subroutines | |
real(kind=f64) function | gauss_legendre_integral_1d (f, a, b, n) |
Gauss-Legendre Quadrature. More... | |
real(kind=f64) function, dimension(2, 1:npoints), public | sll_f_gauss_legendre_points_and_weights (npoints, a, b) |
Returns a 2d array of size (2,npoints) containing gauss-legendre points and weights in the interval [a,b]. More... | |
real(kind=f64) function, dimension(1:npoints), public | sll_f_gauss_legendre_points (npoints, a, b) |
real(kind=f64) function, dimension(1:npoints), public | sll_f_gauss_legendre_weights (npoints, a, b) |
|
private |
Gauss-Legendre Quadrature.
To integrate the function f(x) (real-valued and of a single, real-valued argument x) over the interval [a,b], we use the Gauss-Legendre formula
\[ \int_{-1}^1 f(x)dx \approx \sum_{k=1}^{n} w_k f(x_k) \]
where n represents the desired number of Gauss points.
the function maps the interval [-1,1] into the arbitrary interval [a,b].
To be considered is to split this function into degree-specific functions to avoid the select statement.
f | the function to be integrated | |
[in] | a | left-bound of the definition interval of f |
[in] | b | right-bound of the definition interval of f |
[in] | n | the desired number of Gauss points |
Definition at line 164 of file sll_m_gauss_legendre_integration.F90.
real(kind=f64) function, dimension(1:npoints), public sll_m_gauss_legendre_integration::sll_f_gauss_legendre_points | ( | integer(kind=i32), intent(in) | npoints, |
real(kind=f64), intent(in), optional | a, | ||
real(kind=f64), intent(in), optional | b | ||
) |
Definition at line 300 of file sll_m_gauss_legendre_integration.F90.
real(kind=f64) function, dimension(2, 1:npoints), public sll_m_gauss_legendre_integration::sll_f_gauss_legendre_points_and_weights | ( | integer(kind=i32), intent(in) | npoints, |
real(kind=f64), intent(in), optional | a, | ||
real(kind=f64), intent(in), optional | b | ||
) |
Returns a 2d array of size (2,npoints) containing gauss-legendre points and weights in the interval [a,b].
[in] | npoints | Number of gauss points. |
[in] | a | OPTIONAL Minimum value of the interval. |
[in] | b | OPTIONAL Maximun value of the interval. |
Definition at line 264 of file sll_m_gauss_legendre_integration.F90.
real(kind=f64) function, dimension(1:npoints), public sll_m_gauss_legendre_integration::sll_f_gauss_legendre_weights | ( | integer(kind=i32), intent(in) | npoints, |
real(kind=f64), intent(in), optional | a, | ||
real(kind=f64), intent(in), optional | b | ||
) |
Definition at line 334 of file sll_m_gauss_legendre_integration.F90.