bspline_oo_module Module

Object-oriented style wrappers to bspline_sub_module. This module provides classes (bspline_1d, bspline_2d, bspline_3d, bspline_4d, bspline_5d, and bspline_6d) which can be used instead of the main subroutine interface.


Uses

  • module~~bspline_oo_module~~UsesGraph module~bspline_oo_module bspline_oo_module iso_fortran_env iso_fortran_env module~bspline_oo_module->iso_fortran_env module~bspline_kinds_module bspline_kinds_module module~bspline_oo_module->module~bspline_kinds_module module~bspline_sub_module bspline_sub_module module~bspline_oo_module->module~bspline_sub_module module~bspline_kinds_module->iso_fortran_env module~bspline_sub_module->iso_fortran_env module~bspline_sub_module->module~bspline_kinds_module

Used by

  • module~~bspline_oo_module~~UsedByGraph module~bspline_oo_module bspline_oo_module module~bspline_module bspline_module module~bspline_module->module~bspline_oo_module

Variables

Type Visibility Attributes Name Initial
integer(kind=ip), private, parameter :: int_size = storage_size(1_ip, kind=ip)

size of a default integer [bits]

integer(kind=ip), private, parameter :: logical_size = storage_size(.true., kind=ip)

size of a default logical [bits]

integer(kind=ip), private, parameter :: real_size = storage_size(1.0_wp, kind=ip)

size of a real(wp) [bits]


Interfaces

public interface bspline_1d

Constructor for bspline_1d

  • private pure elemental function bspline_1d_constructor_empty() result(me)

    It returns an empty bspline_1d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_1d)

  • private pure function bspline_1d_constructor_auto_knots(x, fcn, kx, extrap) result(me)

    Constructor for a bspline_1d type (auto knots). This is a wrapper for initialize_1d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: fcn

    (nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_1d)

  • private pure function bspline_1d_constructor_specify_knots(x, fcn, kx, tx, extrap) result(me)

    Constructor for a bspline_1d type (user-specified knots). This is a wrapper for initialize_1d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: fcn

    (nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_1d)

public interface bspline_2d

Constructor for bspline_2d

  • private elemental function bspline_2d_constructor_empty() result(me)

    It returns an empty bspline_2d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_2d)

  • private pure function bspline_2d_constructor_auto_knots(x, y, fcn, kx, ky, extrap) result(me)

    Constructor for a bspline_2d type (auto knots). This is a wrapper for initialize_2d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:) :: fcn

    (nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_2d)

  • private pure function bspline_2d_constructor_specify_knots(x, y, fcn, kx, ky, tx, ty, extrap) result(me)

    Constructor for a bspline_2d type (user-specified knots). This is a wrapper for initialize_2d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:) :: fcn

    (nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ty

    The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_2d)

public interface bspline_3d

Constructor for bspline_3d

  • private elemental function bspline_3d_constructor_empty() result(me)

    It returns an empty bspline_3d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_3d)

  • private pure function bspline_3d_constructor_auto_knots(x, y, z, fcn, kx, ky, kz, extrap) result(me)

    Constructor for a bspline_3d type (auto knots). This is a wrapper for initialize_3d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:) :: fcn

    (nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_3d)

  • private pure function bspline_3d_constructor_specify_knots(x, y, z, fcn, kx, ky, kz, tx, ty, tz, extrap) result(me)

    Constructor for a bspline_3d type (user-specified knots). This is a wrapper for initialize_3d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:) :: fcn

    (nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ty

    The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tz

    The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_3d)

public interface bspline_4d

Constructor for bspline_4d

  • private elemental function bspline_4d_constructor_empty() result(me)

    It returns an empty bspline_4d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_4d)

  • private pure function bspline_4d_constructor_auto_knots(x, y, z, q, fcn, kx, ky, kz, kq, extrap) result(me)

    Constructor for a bspline_4d type (auto knots). This is a wrapper for initialize_4d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

    (nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_4d)

  • private pure function bspline_4d_constructor_specify_knots(x, y, z, q, fcn, kx, ky, kz, kq, tx, ty, tz, tq, extrap) result(me)

    Constructor for a bspline_4d type (user-specified knots). This is a wrapper for initialize_4d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

    (nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ty

    The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tz

    The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tq

    The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_4d)

public interface bspline_5d

Constructor for bspline_5d

  • private elemental function bspline_5d_constructor_empty() result(me)

    It returns an empty bspline_5d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_5d)

  • private pure function bspline_5d_constructor_auto_knots(x, y, z, q, r, fcn, kx, ky, kz, kq, kr, extrap) result(me)

    Constructor for a bspline_5d type (auto knots). This is a wrapper for initialize_5d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: r

    (nr) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

    (nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kr

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_5d)

  • private pure function bspline_5d_constructor_specify_knots(x, y, z, q, r, fcn, kx, ky, kz, kq, kr, tx, ty, tz, tq, tr, extrap) result(me)

    Constructor for a bspline_5d type (user-specified knots). This is a wrapper for initialize_5d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: r

    (nr) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

    (nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kr

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ty

    The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tz

    The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tq

    The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tr

    The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_5d)

public interface bspline_6d

Constructor for bspline_6d

  • private elemental function bspline_6d_constructor_empty() result(me)

    It returns an empty bspline_6d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

    Arguments

    None

    Return Value type(bspline_6d)

  • private pure function bspline_6d_constructor_auto_knots(x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, extrap) result(me)

    Constructor for a bspline_6d type (auto knots). This is a wrapper for initialize_6d_auto_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: r

    (nr) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: s

    (ns) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

    (nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kr

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ks

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_6d)

  • private pure function bspline_6d_constructor_specify_knots(x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, tx, ty, tz, tq, tr, ts, extrap) result(me)

    Constructor for a bspline_6d type (user-specified knots). This is a wrapper for initialize_6d_specify_knots.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=wp), intent(in), dimension(:) :: x

    (nx) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: y

    (ny) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: z

    (nz) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: q

    (nq) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: r

    (nr) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:) :: s

    (ns) array of abcissae. Must be strictly increasing.

    real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

    (nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

    integer(kind=ip), intent(in) :: kx

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ky

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kz

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kq

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: kr

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    integer(kind=ip), intent(in) :: ks

    The order of spline pieces in ( ) (order = polynomial degree + 1)

    real(kind=wp), intent(in), dimension(:) :: tx

    The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ty

    The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tz

    The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tq

    The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: tr

    The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

    real(kind=wp), intent(in), dimension(:) :: ts

    The (ns+ks) knots in the direction for the spline interpolant. Must be non-decreasing.

    logical, intent(in), optional :: extrap

    if true, then extrapolation is allowed (default is false)

    Return Value type(bspline_6d)


Abstract Interfaces

abstract interface

  • private pure function size_func(me) result(s)

    interface for size routines

    Arguments

    Type IntentOptional Attributes Name
    class(bspline_class), intent(in) :: me

    Return Value integer(kind=ip)

    size of the structure in bits

abstract interface

  • private pure subroutine destroy_func(me)

    interface for bspline destructor routines

    Arguments

    Type IntentOptional Attributes Name
    class(bspline_class), intent(inout) :: me

Derived Types

type, public ::  bspline_class

Base class for the b-spline types

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: inbvx = 1_ip

internal variable used by dbvalu for efficient processing

integer(kind=ip), private :: iflag = 1_ip

saved iflag from the list routine call.

logical, private :: initialized = .false.

true if the class is initialized and ready to use

logical, private :: extrap = .false.

if true, then extrapolation is allowed during evaluation

Type-Bound Procedures

procedure, private, non_overridable :: destroy_base ../../

destructor for the abstract type

procedure, private, non_overridable :: set_extrap_flag ../../

internal routine to set the extrap flag

procedure(destroy_func), public, deferred :: destroy ../../

destructor

procedure(size_func), public, deferred :: size_of ../../

size of the structure in bits

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

type, public, extends(bspline_class) ::  bspline_1d

Class for 1d b-spline interpolation.

Read more…

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: work_val_1

[[db1val] work array of dimension 3*kx

Constructor

Constructor for bspline_1d

private pure, elemental function bspline_1d_constructor_empty ()

It returns an empty bspline_1d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_1d_constructor_auto_knots (x, fcn, kx, extrap)

Constructor for a bspline_1d type (auto knots). This is a wrapper for initialize_1d_auto_knots.

private pure function bspline_1d_constructor_specify_knots (x, fcn, kx, tx, extrap)

Constructor for a bspline_1d type (user-specified knots). This is a wrapper for initialize_1d_specify_knots.

Finalizations Procedures

final :: finalize_1d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_1d_auto_knots, initialize_1d_specify_knots
procedure, private :: initialize_1d_auto_knots
procedure, private :: initialize_1d_specify_knots
procedure, public :: evaluate => evaluate_1d
procedure, public :: destroy => destroy_1d
procedure, public :: size_of => size_1d
procedure, public :: integral => integral_1d
procedure, public :: fintegral => fintegral_1d

type, public, extends(bspline_class) ::  bspline_2d

Class for 2d b-spline interpolation.

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: ny = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ky = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:,:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ty

The knots in the direction for the spline interpolant

integer(kind=ip), private :: inbvy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloy = 1_ip

internal variable used for efficient processing

real(kind=wp), private, dimension(:), allocatable :: work_val_1

[[db2val] work array of dimension ky

real(kind=wp), private, dimension(:), allocatable :: work_val_2

[[db2val] work array of dimension 3_ip*max(kx,ky)

Constructor

Constructor for bspline_2d

private elemental function bspline_2d_constructor_empty ()

It returns an empty bspline_2d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_2d_constructor_auto_knots (x, y, fcn, kx, ky, extrap)

Constructor for a bspline_2d type (auto knots). This is a wrapper for initialize_2d_auto_knots.

private pure function bspline_2d_constructor_specify_knots (x, y, fcn, kx, ky, tx, ty, extrap)

Constructor for a bspline_2d type (user-specified knots). This is a wrapper for initialize_2d_specify_knots.

Finalizations Procedures

final :: finalize_2d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_2d_auto_knots, initialize_2d_specify_knots
procedure, private :: initialize_2d_auto_knots
procedure, private :: initialize_2d_specify_knots
procedure, public :: evaluate => evaluate_2d
procedure, public :: destroy => destroy_2d
procedure, public :: size_of => size_2d

type, public, extends(bspline_class) ::  bspline_3d

Class for 3d b-spline interpolation.

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: ny = 0_ip

Number of abcissae

integer(kind=ip), private :: nz = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ky = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kz = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:,:,:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ty

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tz

The knots in the direction for the spline interpolant

integer(kind=ip), private :: inbvy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloz = 1_ip

internal variable used for efficient processing

real(kind=wp), private, dimension(:,:), allocatable :: work_val_1

[[db3val] work array of dimension ky,kz

real(kind=wp), private, dimension(:), allocatable :: work_val_2

[[db3val] work array of dimension kz

real(kind=wp), private, dimension(:), allocatable :: work_val_3

[[db3val] work array of dimension 3_ip*max(kx,ky,kz)

Constructor

Constructor for bspline_3d

private elemental function bspline_3d_constructor_empty ()

It returns an empty bspline_3d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_3d_constructor_auto_knots (x, y, z, fcn, kx, ky, kz, extrap)

Constructor for a bspline_3d type (auto knots). This is a wrapper for initialize_3d_auto_knots.

private pure function bspline_3d_constructor_specify_knots (x, y, z, fcn, kx, ky, kz, tx, ty, tz, extrap)

Constructor for a bspline_3d type (user-specified knots). This is a wrapper for initialize_3d_specify_knots.

Finalizations Procedures

final :: finalize_3d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_3d_auto_knots, initialize_3d_specify_knots
procedure, private :: initialize_3d_auto_knots
procedure, private :: initialize_3d_specify_knots
procedure, public :: evaluate => evaluate_3d
procedure, public :: destroy => destroy_3d
procedure, public :: size_of => size_3d

type, public, extends(bspline_class) ::  bspline_4d

Class for 4d b-spline interpolation.

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: ny = 0_ip

Number of abcissae

integer(kind=ip), private :: nz = 0_ip

Number of abcissae

integer(kind=ip), private :: nq = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ky = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kz = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kq = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:,:,:,:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ty

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tz

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tq

The knots in the direction for the spline interpolant

integer(kind=ip), private :: inbvy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvq = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloq = 1_ip

internal variable used for efficient processing

real(kind=wp), private, dimension(:,:,:), allocatable :: work_val_1

db4val work array of dimension ky,kz,kq

real(kind=wp), private, dimension(:,:), allocatable :: work_val_2

db4val work array of dimension kz,kq

real(kind=wp), private, dimension(:), allocatable :: work_val_3

db4val work array of dimension kq

real(kind=wp), private, dimension(:), allocatable :: work_val_4

db4val work array of dimension 3_ip*max(kx,ky,kz,kq)

Constructor

Constructor for bspline_4d

private elemental function bspline_4d_constructor_empty ()

It returns an empty bspline_4d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_4d_constructor_auto_knots (x, y, z, q, fcn, kx, ky, kz, kq, extrap)

Constructor for a bspline_4d type (auto knots). This is a wrapper for initialize_4d_auto_knots.

private pure function bspline_4d_constructor_specify_knots (x, y, z, q, fcn, kx, ky, kz, kq, tx, ty, tz, tq, extrap)

Constructor for a bspline_4d type (user-specified knots). This is a wrapper for initialize_4d_specify_knots.

Finalizations Procedures

final :: finalize_4d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_4d_auto_knots, initialize_4d_specify_knots
procedure, private :: initialize_4d_auto_knots
procedure, private :: initialize_4d_specify_knots
procedure, public :: evaluate => evaluate_4d
procedure, public :: destroy => destroy_4d
procedure, public :: size_of => size_4d

type, public, extends(bspline_class) ::  bspline_5d

Class for 5d b-spline interpolation.

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: ny = 0_ip

Number of abcissae

integer(kind=ip), private :: nz = 0_ip

Number of abcissae

integer(kind=ip), private :: nq = 0_ip

Number of abcissae

integer(kind=ip), private :: nr = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ky = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kz = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kq = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kr = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:,:,:,:,:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ty

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tz

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tq

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tr

The knots in the direction for the spline interpolant

integer(kind=ip), private :: inbvy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvq = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvr = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloq = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: ilor = 1_ip

internal variable used for efficient processing

real(kind=wp), private, dimension(:,:,:,:), allocatable :: work_val_1

db5val work array of dimension ky,kz,kq,kr

real(kind=wp), private, dimension(:,:,:), allocatable :: work_val_2

db5val work array of dimension kz,kq,kr

real(kind=wp), private, dimension(:,:), allocatable :: work_val_3

db5val work array of dimension kq,kr

real(kind=wp), private, dimension(:), allocatable :: work_val_4

db5val work array of dimension kr

real(kind=wp), private, dimension(:), allocatable :: work_val_5

db5val work array of dimension 3_ip*max(kx,ky,kz,kq,kr)

Constructor

Constructor for bspline_5d

private elemental function bspline_5d_constructor_empty ()

It returns an empty bspline_5d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_5d_constructor_auto_knots (x, y, z, q, r, fcn, kx, ky, kz, kq, kr, extrap)

Constructor for a bspline_5d type (auto knots). This is a wrapper for initialize_5d_auto_knots.

private pure function bspline_5d_constructor_specify_knots (x, y, z, q, r, fcn, kx, ky, kz, kq, kr, tx, ty, tz, tq, tr, extrap)

Constructor for a bspline_5d type (user-specified knots). This is a wrapper for initialize_5d_specify_knots.

Finalizations Procedures

final :: finalize_5d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_5d_auto_knots, initialize_5d_specify_knots
procedure, private :: initialize_5d_auto_knots
procedure, private :: initialize_5d_specify_knots
procedure, public :: evaluate => evaluate_5d
procedure, public :: destroy => destroy_5d
procedure, public :: size_of => size_5d

type, public, extends(bspline_class) ::  bspline_6d

Class for 6d b-spline interpolation.

Components

Type Visibility Attributes Name Initial
integer(kind=ip), private :: nx = 0_ip

Number of abcissae

integer(kind=ip), private :: ny = 0_ip

Number of abcissae

integer(kind=ip), private :: nz = 0_ip

Number of abcissae

integer(kind=ip), private :: nq = 0_ip

Number of abcissae

integer(kind=ip), private :: nr = 0_ip

Number of abcissae

integer(kind=ip), private :: ns = 0_ip

Number of abcissae

integer(kind=ip), private :: kx = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ky = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kz = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kq = 0_ip

The order of spline pieces in

integer(kind=ip), private :: kr = 0_ip

The order of spline pieces in

integer(kind=ip), private :: ks = 0_ip

The order of spline pieces in

real(kind=wp), private, dimension(:,:,:,:,:,:), allocatable :: bcoef

array of coefficients of the b-spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tx

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ty

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tz

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tq

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: tr

The knots in the direction for the spline interpolant

real(kind=wp), private, dimension(:), allocatable :: ts

The knots in the direction for the spline interpolant

integer(kind=ip), private :: inbvy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvq = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvr = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: inbvs = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloy = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloz = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: iloq = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: ilor = 1_ip

internal variable used for efficient processing

integer(kind=ip), private :: ilos = 1_ip

internal variable used for efficient processing

real(kind=wp), private, dimension(:,:,:,:,:), allocatable :: work_val_1

db6val work array of dimension ky,kz,kq,kr,ks

real(kind=wp), private, dimension(:,:,:,:), allocatable :: work_val_2

db6val work array of dimension kz,kq,kr,ks

real(kind=wp), private, dimension(:,:,:), allocatable :: work_val_3

db6val work array of dimension kq,kr,ks

real(kind=wp), private, dimension(:,:), allocatable :: work_val_4

db6val work array of dimension kr,ks

real(kind=wp), private, dimension(:), allocatable :: work_val_5

db6val work array of dimension ks

real(kind=wp), private, dimension(:), allocatable :: work_val_6

db6val work array of dimension 3_ip*max(kx,ky,kz,kq,kr,ks)

Constructor

Constructor for bspline_6d

private elemental function bspline_6d_constructor_empty ()

It returns an empty bspline_6d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

private pure function bspline_6d_constructor_auto_knots (x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, extrap)

Constructor for a bspline_6d type (auto knots). This is a wrapper for initialize_6d_auto_knots.

private pure function bspline_6d_constructor_specify_knots (x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, tx, ty, tz, tq, tr, ts, extrap)

Constructor for a bspline_6d type (user-specified knots). This is a wrapper for initialize_6d_specify_knots.

Finalizations Procedures

final :: finalize_6d

Type-Bound Procedures

procedure, public, non_overridable :: status_ok ../../

returns true if the last iflag status code was =0.

procedure, public, non_overridable :: status_message => get_bspline_status_message ../../

retrieve the last status message

procedure, public, non_overridable :: clear_flag => clear_bspline_flag ../../

to reset the iflag saved in the class.

generic, public :: initialize => initialize_6d_auto_knots, initialize_6d_specify_knots
procedure, private :: initialize_6d_auto_knots
procedure, private :: initialize_6d_specify_knots
procedure, public :: evaluate => evaluate_6d
procedure, public :: destroy => destroy_6d
procedure, public :: size_of => size_6d

Functions

private elemental function status_ok(me) result(ok)

This routines returns true if the iflag code from the last routine called was =0. Maybe of the routines have output iflag variables, so they can be checked explicitly, or this routine can be used.

Read more…

Arguments

Type IntentOptional Attributes Name
class(bspline_class), intent(in) :: me

Return Value logical

private pure function get_bspline_status_message(me, iflag) result(msg)

Get the status message from a bspline_class routine call.

Read more…

Arguments

Type IntentOptional Attributes Name
class(bspline_class), intent(in) :: me
integer(kind=ip), intent(in), optional :: iflag

the corresponding status code

Return Value character(len=:), allocatable

status message associated with the flag

private pure function size_1d(me) result(s)

Actual size of a bspline_1d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure function size_2d(me) result(s)

Actual size of a bspline_2d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_2d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure function size_3d(me) result(s)

Actual size of a bspline_3d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_3d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure function size_4d(me) result(s)

Actual size of a bspline_4d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_4d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure function size_5d(me) result(s)

Actual size of a bspline_5d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_5d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure function size_6d(me) result(s)

Actual size of a bspline_6d structure in bits.

Arguments

Type IntentOptional Attributes Name
class(bspline_6d), intent(in) :: me

Return Value integer(kind=ip)

size of the structure in bits

private pure elemental function bspline_1d_constructor_empty() result(me)

It returns an empty bspline_1d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_1d)

private pure function bspline_1d_constructor_auto_knots(x, fcn, kx, extrap) result(me)

Constructor for a bspline_1d type (auto knots). This is a wrapper for initialize_1d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: fcn

(nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_1d)

private pure function bspline_1d_constructor_specify_knots(x, fcn, kx, tx, extrap) result(me)

Constructor for a bspline_1d type (user-specified knots). This is a wrapper for initialize_1d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: fcn

(nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_1d)

private elemental function bspline_2d_constructor_empty() result(me)

It returns an empty bspline_2d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_2d)

private pure function bspline_2d_constructor_auto_knots(x, y, fcn, kx, ky, extrap) result(me)

Constructor for a bspline_2d type (auto knots). This is a wrapper for initialize_2d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:) :: fcn

(nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_2d)

private pure function bspline_2d_constructor_specify_knots(x, y, fcn, kx, ky, tx, ty, extrap) result(me)

Constructor for a bspline_2d type (user-specified knots). This is a wrapper for initialize_2d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:) :: fcn

(nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_2d)

private elemental function bspline_3d_constructor_empty() result(me)

It returns an empty bspline_3d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_3d)

private pure function bspline_3d_constructor_auto_knots(x, y, z, fcn, kx, ky, kz, extrap) result(me)

Constructor for a bspline_3d type (auto knots). This is a wrapper for initialize_3d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:) :: fcn

(nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_3d)

private pure function bspline_3d_constructor_specify_knots(x, y, z, fcn, kx, ky, kz, tx, ty, tz, extrap) result(me)

Constructor for a bspline_3d type (user-specified knots). This is a wrapper for initialize_3d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:) :: fcn

(nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_3d)

private elemental function bspline_4d_constructor_empty() result(me)

It returns an empty bspline_4d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_4d)

private pure function bspline_4d_constructor_auto_knots(x, y, z, q, fcn, kx, ky, kz, kq, extrap) result(me)

Constructor for a bspline_4d type (auto knots). This is a wrapper for initialize_4d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

(nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_4d)

private pure function bspline_4d_constructor_specify_knots(x, y, z, q, fcn, kx, ky, kz, kq, tx, ty, tz, tq, extrap) result(me)

Constructor for a bspline_4d type (user-specified knots). This is a wrapper for initialize_4d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

(nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_4d)

private elemental function bspline_5d_constructor_empty() result(me)

It returns an empty bspline_5d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_5d)

private pure function bspline_5d_constructor_auto_knots(x, y, z, q, r, fcn, kx, ky, kz, kq, kr, extrap) result(me)

Constructor for a bspline_5d type (auto knots). This is a wrapper for initialize_5d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_5d)

private pure function bspline_5d_constructor_specify_knots(x, y, z, q, r, fcn, kx, ky, kz, kq, kr, tx, ty, tz, tq, tr, extrap) result(me)

Constructor for a bspline_5d type (user-specified knots). This is a wrapper for initialize_5d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tr

The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_5d)

private elemental function bspline_6d_constructor_empty() result(me)

It returns an empty bspline_6d type. Note that INITIALIZE still needs to be called before it can be used. Not really that useful except perhaps in some OpenMP applications.

Arguments

None

Return Value type(bspline_6d)

private pure function bspline_6d_constructor_auto_knots(x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, extrap) result(me)

Constructor for a bspline_6d type (auto knots). This is a wrapper for initialize_6d_auto_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: s

(ns) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ks

The order of spline pieces in ( ) (order = polynomial degree + 1)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_6d)

private pure function bspline_6d_constructor_specify_knots(x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, tx, ty, tz, tq, tr, ts, extrap) result(me)

Constructor for a bspline_6d type (user-specified knots). This is a wrapper for initialize_6d_specify_knots.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: s

(ns) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ks

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tr

The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ts

The (ns+ks) knots in the direction for the spline interpolant. Must be non-decreasing.

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

Return Value type(bspline_6d)


Subroutines

private elemental subroutine clear_bspline_flag(me)

This sets the iflag variable in the class to 0 (which indicates that everything is OK). It can be used after an error is encountered.

Arguments

Type IntentOptional Attributes Name
class(bspline_class), intent(inout) :: me

private pure subroutine destroy_base(me)

Destructor for contents of the base bspline_class class. (this routine is called by the extended classes).

Arguments

Type IntentOptional Attributes Name
class(bspline_class), intent(inout) :: me

private pure subroutine destroy_1d(me)

Destructor for bspline_1d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me

private pure subroutine destroy_2d(me)

Destructor for bspline_2d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_2d), intent(inout) :: me

private pure subroutine destroy_3d(me)

Destructor for bspline_3d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_3d), intent(inout) :: me

private pure subroutine destroy_4d(me)

Destructor for bspline_4d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_4d), intent(inout) :: me

private pure subroutine destroy_5d(me)

Destructor for bspline_5d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_5d), intent(inout) :: me

private pure subroutine destroy_6d(me)

Destructor for bspline_6d class.

Arguments

Type IntentOptional Attributes Name
class(bspline_6d), intent(inout) :: me

private pure elemental subroutine finalize_1d(me)

Finalizer for bspline_1d class. Just a wrapper for destroy_1d.

Arguments

Type IntentOptional Attributes Name
type(bspline_1d), intent(inout) :: me

private pure elemental subroutine finalize_2d(me)

Finalizer for bspline_2d class. Just a wrapper for destroy_2d.

Arguments

Type IntentOptional Attributes Name
type(bspline_2d), intent(inout) :: me

private pure elemental subroutine finalize_3d(me)

Finalizer for bspline_3d class. Just a wrapper for destroy_3d.

Arguments

Type IntentOptional Attributes Name
type(bspline_3d), intent(inout) :: me

private pure elemental subroutine finalize_4d(me)

Finalizer for bspline_4d class. Just a wrapper for destroy_4d.

Arguments

Type IntentOptional Attributes Name
type(bspline_4d), intent(inout) :: me

private pure elemental subroutine finalize_5d(me)

Finalizer for bspline_5d class. Just a wrapper for destroy_5d.

Arguments

Type IntentOptional Attributes Name
type(bspline_5d), intent(inout) :: me

private pure elemental subroutine finalize_6d(me)

Finalizer for bspline_6d class. Just a wrapper for destroy_6d.

Arguments

Type IntentOptional Attributes Name
type(bspline_6d), intent(inout) :: me

private pure subroutine set_extrap_flag(me, extrap)

Sets the extrap flag in the class.

Arguments

Type IntentOptional Attributes Name
class(bspline_class), intent(inout) :: me
logical, intent(in), optional :: extrap

if not present, then False is used

private pure subroutine initialize_1d_auto_knots(me, x, fcn, kx, iflag, extrap)

Initialize a bspline_1d type (with automatically-computed knots). This is a wrapper for db1ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: fcn

(nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db1ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_1d_specify_knots(me, x, fcn, kx, tx, iflag, extrap)

Initialize a bspline_1d type (with user-specified knots). This is a wrapper for db1ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: fcn

(nx) array of function values to interpolate. fcn(i) should contain the function value at the point x(i)

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db1ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_1d(me, xval, idx, f, iflag)

Evaluate a bspline_1d interpolate. This is a wrapper for db1val.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db1val)

private pure subroutine integral_1d(me, x1, x2, f, iflag)

Evaluate a bspline_1d definite integral. This is a wrapper for db1sqad.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me
real(kind=wp), intent(in) :: x1

left point of interval

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

right point of interval

real(kind=wp), intent(out) :: f

integral of the b-spline over

integer(kind=ip), intent(out) :: iflag

status flag (see db1sqad)

private subroutine fintegral_1d(me, fun, idx, x1, x2, tol, f, iflag)

Evaluate a bspline_1d definite integral. This is a wrapper for db1fqad.

Arguments

Type IntentOptional Attributes Name
class(bspline_1d), intent(inout) :: me
procedure(b1fqad_func) :: fun

external function of one argument for the integrand bf(x)=fun(x)*dbvalu(tx,bcoef,nx,kx,idx,x,inbv)

integer(kind=ip), intent(in) :: idx

order of the spline derivative, 0 <= idx <= k-1 idx=0 gives the spline function

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

left point of interval

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

right point of interval

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

desired accuracy for the quadrature

real(kind=wp), intent(out) :: f

integral of bf(x) over

integer(kind=ip), intent(out) :: iflag

status flag (see db1sqad)

private pure subroutine initialize_2d_auto_knots(me, x, y, fcn, kx, ky, iflag, extrap)

Initialize a bspline_2d type (with automatically-computed knots). This is a wrapper for db2ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_2d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:) :: fcn

(nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db2ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_2d_specify_knots(me, x, y, fcn, kx, ky, tx, ty, iflag, extrap)

Initialize a bspline_2d type (with user-specified knots). This is a wrapper for db2ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_2d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:) :: fcn

(nx,ny) matrix of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db2ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_2d(me, xval, yval, idx, idy, f, iflag)

Evaluate a bspline_2d interpolate. This is a wrapper for db2val.

Arguments

Type IntentOptional Attributes Name
class(bspline_2d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

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

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idy

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db2val)

private pure subroutine initialize_3d_auto_knots(me, x, y, z, fcn, kx, ky, kz, iflag, extrap)

Initialize a bspline_3d type (with automatically-computed knots). This is a wrapper for db3ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_3d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:) :: fcn

(nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db3ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_3d_specify_knots(me, x, y, z, fcn, kx, ky, kz, tx, ty, tz, iflag, extrap)

Initialize a bspline_3d type (with user-specified knots). This is a wrapper for db3ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_3d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:) :: fcn

(nx,ny,nz) matrix of function values to interpolate. fcn(i,j,k) should contain the function value at the point (x(i),y(j),z(k))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db3ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_3d(me, xval, yval, zval, idx, idy, idz, f, iflag)

Evaluate a bspline_3d interpolate. This is a wrapper for db3val.

Arguments

Type IntentOptional Attributes Name
class(bspline_3d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idy

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idz

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db3val)

private pure subroutine initialize_4d_auto_knots(me, x, y, z, q, fcn, kx, ky, kz, kq, iflag, extrap)

Initialize a bspline_4d type (with automatically-computed knots). This is a wrapper for db4ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_4d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

(nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db4ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_4d_specify_knots(me, x, y, z, q, fcn, kx, ky, kz, kq, tx, ty, tz, tq, iflag, extrap)

Initialize a bspline_4d type (with user-specified knots). This is a wrapper for db4ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_4d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:) :: fcn

(nx,ny,nz,nq) matrix of function values to interpolate. fcn(i,j,k,l) should contain the function value at the point (x(i),y(j),z(k),q(l))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db4ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_4d(me, xval, yval, zval, qval, idx, idy, idz, idq, f, iflag)

Evaluate a bspline_4d interpolate. This is a wrapper for db4val.

Arguments

Type IntentOptional Attributes Name
class(bspline_4d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idy

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idz

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idq

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db4val)

private pure subroutine initialize_5d_auto_knots(me, x, y, z, q, r, fcn, kx, ky, kz, kq, kr, iflag, extrap)

Initialize a bspline_5d type (with automatically-computed knots). This is a wrapper for db5ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_5d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db5ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_5d_specify_knots(me, x, y, z, q, r, fcn, kx, ky, kz, kq, kr, tx, ty, tz, tq, tr, iflag, extrap)

Initialize a bspline_5d type (with user-specified knots). This is a wrapper for db5ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_5d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr) matrix of function values to interpolate. fcn(i,j,k,l,m) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tr

The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db5ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_5d(me, xval, yval, zval, qval, rval, idx, idy, idz, idq, idr, f, iflag)

Evaluate a bspline_5d interpolate. This is a wrapper for db5val.

Arguments

Type IntentOptional Attributes Name
class(bspline_5d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idy

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idz

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idq

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idr

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db5val)

private pure subroutine initialize_6d_auto_knots(me, x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, iflag, extrap)

Initialize a bspline_6d type (with automatically-computed knots). This is a wrapper for db6ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_6d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: s

(ns) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ks

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(out) :: iflag

status flag (see db6ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine initialize_6d_specify_knots(me, x, y, z, q, r, s, fcn, kx, ky, kz, kq, kr, ks, tx, ty, tz, tq, tr, ts, iflag, extrap)

Initialize a bspline_6d type (with user-specified knots). This is a wrapper for db6ink.

Arguments

Type IntentOptional Attributes Name
class(bspline_6d), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

(nx) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: y

(ny) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: z

(nz) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: q

(nq) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: r

(nr) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:) :: s

(ns) array of abcissae. Must be strictly increasing.

real(kind=wp), intent(in), dimension(:,:,:,:,:,:) :: fcn

(nx,ny,nz,nq,nr,ns) matrix of function values to interpolate. fcn(i,j,k,l,m,n) should contain the function value at the point (x(i),y(j),z(k),q(l),r(m),s(n))

integer(kind=ip), intent(in) :: kx

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ky

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kz

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kq

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: kr

The order of spline pieces in ( ) (order = polynomial degree + 1)

integer(kind=ip), intent(in) :: ks

The order of spline pieces in ( ) (order = polynomial degree + 1)

real(kind=wp), intent(in), dimension(:) :: tx

The (nx+kx) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ty

The (ny+ky) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tz

The (nz+kz) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tq

The (nq+kq) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: tr

The (nr+kr) knots in the direction for the spline interpolant. Must be non-decreasing.

real(kind=wp), intent(in), dimension(:) :: ts

The (ns+ks) knots in the direction for the spline interpolant. Must be non-decreasing.

integer(kind=ip), intent(out) :: iflag

status flag (see db6ink)

logical, intent(in), optional :: extrap

if true, then extrapolation is allowed (default is false)

private pure subroutine evaluate_6d(me, xval, yval, zval, qval, rval, sval, idx, idy, idz, idq, idr, ids, f, iflag)

Evaluate a bspline_6d interpolate. This is a wrapper for db6val.

Arguments

Type IntentOptional Attributes Name
class(bspline_6d), intent(inout) :: me
real(kind=wp), intent(in) :: xval

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

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

coordinate of evaluation point.

integer(kind=ip), intent(in) :: idx

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idy

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idz

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idq

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: idr

derivative of piecewise polynomial to evaluate.

integer(kind=ip), intent(in) :: ids

derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(out) :: f

interpolated value

integer(kind=ip), intent(out) :: iflag

status flag (see db6val)

private pure subroutine check_knot_vectors_sizes(nx, ny, nz, nq, nr, ns, kx, ky, kz, kq, kr, ks, tx, ty, tz, tq, tr, ts, iflag)

Error checks for the user-specified knot vector sizes.

Read more…

Arguments

Type IntentOptional Attributes Name
integer(kind=ip), intent(in), optional :: nx
integer(kind=ip), intent(in), optional :: ny
integer(kind=ip), intent(in), optional :: nz
integer(kind=ip), intent(in), optional :: nq
integer(kind=ip), intent(in), optional :: nr
integer(kind=ip), intent(in), optional :: ns
integer(kind=ip), intent(in), optional :: kx
integer(kind=ip), intent(in), optional :: ky
integer(kind=ip), intent(in), optional :: kz
integer(kind=ip), intent(in), optional :: kq
integer(kind=ip), intent(in), optional :: kr
integer(kind=ip), intent(in), optional :: ks
real(kind=wp), intent(in), optional, dimension(:) :: tx
real(kind=wp), intent(in), optional, dimension(:) :: ty
real(kind=wp), intent(in), optional, dimension(:) :: tz
real(kind=wp), intent(in), optional, dimension(:) :: tq
real(kind=wp), intent(in), optional, dimension(:) :: tr
real(kind=wp), intent(in), optional, dimension(:) :: ts
integer(kind=ip), intent(out) :: iflag

0 if everything is OK