Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_rectangle_integration Module Reference

Description

Rectangle integration.

Low-level mathematical utility that applies the Rectangle method to compute numeric integrals.

Derived types and interfaces

interface  sll_o_rectangle_integrate_1d
 Integrate numerically with Gauss-Lobatto formula. More...
 

Functions/Subroutines

real(kind=f64) function rectangle_integral_1d (f, x, n)
 Integrate with rectangle formula. More...
 
real(kind=f64) function, dimension(n) rectangle_weights (n, x)
 Returns a 1d array of size (n) containing rectangle integration weights in the interval [x(1),x(n)]. More...
 

Function/Subroutine Documentation

◆ rectangle_integral_1d()

real(kind=f64) function sll_m_rectangle_integration::rectangle_integral_1d ( procedure(function_1d)  f,
real(kind=f64), dimension(n)  x,
integer(kind=i32), intent(in)  n 
)
private

Integrate with rectangle formula.

To integrate the function \( f(x) \) , we use the rectangle formula

\[ \int_{-1}^1 f(x)dx \approx \sum_{k=1}^{n} w_k f(x_k) \]

where n and x represents the desired number of points and their positions.

Parameters
fthe function to be integrated
[in]xpositions of points
[in]nthe number of points
Returns
The value of the integral

Definition at line 51 of file sll_m_rectangle_integration.F90.

Here is the caller graph for this function:

◆ rectangle_weights()

real(kind=f64) function, dimension(n) sll_m_rectangle_integration::rectangle_weights ( integer(kind=i32), intent(in)  n,
real(kind=f64), dimension(n)  x 
)
private

Returns a 1d array of size (n) containing rectangle integration weights in the interval [x(1),x(n)].

Parameters
[in]nNumber of gauss points.
[in]xPoint poisitions in interval.
Returns
w Array containing weights.

Definition at line 72 of file sll_m_rectangle_integration.F90.

    Report Typos and Errors