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

Description

Access point to B-splines of arbitrary degree providing factory function.

Author
Yaman Güçlü - IPP Garching

Functions/Subroutines

subroutine, public sll_s_bsplines_new (bsplines, degree, periodic, xmin, xmax, ncells, breaks)
 Allocate and initialize uniform or non-uniform B-splines. More...
 
subroutine, public sll_s_bsplines_new_2d_polar (bsplines_radial, bsplines_angular, degree, ncells, max_radius, theta_lims, breaks_radial, breaks_angular)
 Allocate and initialize two B-splines for use on polar grid. More...
 
subroutine, public sll_s_bsplines_new_mirror_copy (radial_basis, extended_basis)
 Create B-splines by mirroring radial basis onto extended domain [-R,R]. More...
 

Variables

integer, parameter wp = f64
 Working precision. More...
 

Function/Subroutine Documentation

◆ sll_s_bsplines_new()

subroutine, public sll_m_bsplines::sll_s_bsplines_new ( class(sll_c_bsplines), intent(inout), allocatable  bsplines,
integer, intent(in)  degree,
logical, intent(in)  periodic,
real(wp), intent(in)  xmin,
real(wp), intent(in)  xmax,
integer, intent(in)  ncells,
real(wp), dimension(:), intent(in), optional  breaks 
)

Allocate and initialize uniform or non-uniform B-splines.

Parameters
[out]bsplinesallocatable polymorphic object
[in]degreespline degree
[in]periodic.true. if domain is periodic
[in]xminx coordinate of left boundary of domain
[in]xmaxx coordinate of right boundary of domain
[in]ncellsnumber of cells in domain (one polynomial per cell)
[in]breakslist of breakpoints (only for non-uniform B-splines)

Definition at line 43 of file sll_m_bsplines.F90.

Here is the caller graph for this function:

◆ sll_s_bsplines_new_2d_polar()

subroutine, public sll_m_bsplines::sll_s_bsplines_new_2d_polar ( class(sll_c_bsplines), intent(inout), allocatable  bsplines_radial,
class(sll_c_bsplines), intent(inout), allocatable  bsplines_angular,
integer, dimension(2), intent(in)  degree,
integer, dimension(2), intent(in)  ncells,
real(wp), intent(in)  max_radius,
real(wp), dimension(2), intent(in)  theta_lims,
real(wp), dimension(:), intent(in), optional  breaks_radial,
real(wp), dimension(:), intent(in), optional  breaks_angular 
)

Allocate and initialize two B-splines for use on polar grid.

Definition at line 101 of file sll_m_bsplines.F90.

Here is the call graph for this function:

◆ sll_s_bsplines_new_mirror_copy()

subroutine, public sll_m_bsplines::sll_s_bsplines_new_mirror_copy ( class(sll_c_bsplines), intent(in)  radial_basis,
class(sll_c_bsplines), intent(inout), allocatable  extended_basis 
)

Create B-splines by mirroring radial basis onto extended domain [-R,R].

Definition at line 181 of file sll_m_bsplines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ wp

integer, parameter wp = f64
private

Working precision.

Definition at line 27 of file sll_m_bsplines.F90.

    Report Typos and Errors