sll_c_bsplines Derived Type

type, public, abstract :: sll_c_bsplines

Abstract type, B-splines


Inherited by

type~~sll_c_bsplines~~InheritedByGraph type~sll_c_bsplines sll_c_bsplines type~sll_t_bsplines_non_uniform sll_t_bsplines_non_uniform type~sll_t_bsplines_non_uniform->type~sll_c_bsplines type~sll_t_bsplines_uniform sll_t_bsplines_uniform type~sll_t_bsplines_uniform->type~sll_c_bsplines type~sll_t_spline_1d sll_t_spline_1d type~sll_t_spline_1d->type~sll_c_bsplines bspl type~sll_t_spline_interpolator_1d sll_t_spline_interpolator_1d type~sll_t_spline_interpolator_1d->type~sll_c_bsplines bspl type~t_cubspl t_cubspl type~t_cubspl->type~sll_c_bsplines bspl type~t_rprofile_bspl t_rProfile_bspl type~t_rprofile_bspl->type~sll_c_bsplines bspl type~t_sbase_spl t_sBase_spl type~t_sbase_spl->type~sll_c_bsplines bspl type~t_sbase_spl->type~sll_t_spline_interpolator_1d Interpol

Components

Type Visibility Attributes Name Initial
integer, public :: degree
logical, public :: periodic
logical, public :: uniform
integer, public :: ncells
integer, public :: nbasis
integer, public :: offset
real(kind=wp), public :: xmin
real(kind=wp), public :: xmax
real(kind=wp), public, allocatable :: knots(:)

Type-Bound Procedures

procedure(i_fun_find_cell), public, deferred :: find_cell

  • pure function i_fun_find_cell(self, x) result(icell) Prototype

    results cell index

    Arguments

    Type IntentOptional Attributes Name
    class(sll_c_bsplines), intent(in) :: self

    B-splines object

    real(kind=wp), intent(in) :: x

    point of interest

    Return Value integer

procedure(i_sub_eval_basis), public, deferred :: eval_basis

  • pure subroutine i_sub_eval_basis(self, x, values, jmin) Prototype

    Evaluate value at x of all basis functions with support in local cell values[j] = B_j(x) for jmin <= j <= jmin+degree

    Arguments

    Type IntentOptional Attributes Name
    class(sll_c_bsplines), intent(in) :: self

    B-splines object

    real(kind=wp), intent(in) :: x

    evaluation point

    real(kind=wp), intent(out) :: values(:)
    integer, intent(out) :: jmin

    index of first non-zero B-spline

procedure(i_sub_eval_deriv), public, deferred :: eval_deriv

  • pure subroutine i_sub_eval_deriv(self, x, derivs, jmin) Prototype

    Evaluate derivative at x of all basis functions with support in local cell derivs[j] = B_j'(x) for jmin <= j <= jmin+degree

    Arguments

    Type IntentOptional Attributes Name
    class(sll_c_bsplines), intent(in) :: self

    B-splines object

    real(kind=wp), intent(in) :: x

    evaluation point

    real(kind=wp), intent(out) :: derivs(:)
    integer, intent(out) :: jmin

    index of first non-zero B-spline

procedure(i_sub_eval_basis_and_n_derivs), public, deferred :: eval_basis_and_n_derivs

  • pure subroutine i_sub_eval_basis_and_n_derivs(self, x, n, derivs, jmin) Prototype

    Evaluate value and n derivatives at x of all basis functions with support in local cell derivs[i,j] = (d/dx)^i B_j(x) for 0 <= i <= n and jmin <= j <= jmin+degree

    Arguments

    Type IntentOptional Attributes Name
    class(sll_c_bsplines), intent(in) :: self

    B-splines object

    real(kind=wp), intent(in) :: x

    evaluation point

    integer, intent(in) :: n

    number of required derivatives

    real(kind=wp), intent(out) :: derivs(:,:)
    integer, intent(out) :: jmin

    index of first non-zero B-spline

procedure(i_sub_free), public, deferred :: free

  • subroutine i_sub_free(self) Prototype

    Arguments

    Type IntentOptional Attributes Name
    class(sll_c_bsplines), intent(inout) :: self

    B-splines object

Source Code

  type, abstract :: sll_c_bsplines

    integer :: degree
    logical :: periodic
    logical :: uniform
    integer :: ncells
    integer :: nbasis
    integer :: offset

    real(wp) :: xmin
    real(wp) :: xmax

    real(wp), allocatable :: knots(:) ! Only used by non-uniform B-splines

  contains
    procedure(i_fun_find_cell              ), deferred :: find_cell
    procedure(i_sub_eval_basis             ), deferred :: eval_basis
    procedure(i_sub_eval_deriv             ), deferred :: eval_deriv
    procedure(i_sub_eval_basis_and_n_derivs), deferred :: eval_basis_and_n_derivs
    procedure(i_sub_free                   ), deferred :: free

  end type sll_c_bsplines