Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Access point to B-splines of arbitrary degree providing factory function.
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... | |
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.
[out] | bsplines | allocatable polymorphic object |
[in] | degree | spline degree |
[in] | periodic | .true. if domain is periodic |
[in] | xmin | x coordinate of left boundary of domain |
[in] | xmax | x coordinate of right boundary of domain |
[in] | ncells | number of cells in domain (one polynomial per cell) |
[in] | breaks | list of breakpoints (only for non-uniform B-splines) |
Definition at line 43 of file sll_m_bsplines.F90.
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.
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.
|
private |
Working precision.
Definition at line 27 of file sll_m_bsplines.F90.