bspline_6d Derived Type

type, public, extends(bspline_class) :: bspline_6d

Class for 6d b-spline interpolation.


Inherits

type~~bspline_6d~~InheritsGraph type~bspline_6d bspline_6d type~bspline_class bspline_class type~bspline_6d->type~bspline_class

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

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)


Finalization Procedures

final :: finalize_6d


Type-Bound Procedures

procedure, public, non_overridable :: status_ok

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

  • 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

procedure, public, non_overridable :: status_message => get_bspline_status_message

retrieve the last status message

  • 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

procedure, public, non_overridable :: clear_flag => clear_bspline_flag

to reset the iflag saved in the class.

  • 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

generic, public :: initialize => initialize_6d_auto_knots, initialize_6d_specify_knots

  • 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)

procedure, private :: initialize_6d_auto_knots

  • 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)

procedure, private :: initialize_6d_specify_knots

  • 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)

procedure, public :: evaluate => evaluate_6d

  • 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)

procedure, public :: destroy => destroy_6d

  • private pure subroutine destroy_6d(me)

    Destructor for bspline_6d class.

    Arguments

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

procedure, public :: size_of => size_6d

  • 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

Source Code

    type,extends(bspline_class),public :: bspline_6d
        !! Class for 6d b-spline interpolation.
        private
        integer(ip) :: nx = 0_ip  !! Number of \(x\) abcissae
        integer(ip) :: ny = 0_ip  !! Number of \(y\) abcissae
        integer(ip) :: nz = 0_ip  !! Number of \(z\) abcissae
        integer(ip) :: nq = 0_ip  !! Number of \(q\) abcissae
        integer(ip) :: nr = 0_ip  !! Number of \(r\) abcissae
        integer(ip) :: ns = 0_ip  !! Number of \(s\) abcissae
        integer(ip) :: kx = 0_ip  !! The order of spline pieces in \(x\)
        integer(ip) :: ky = 0_ip  !! The order of spline pieces in \(y\)
        integer(ip) :: kz = 0_ip  !! The order of spline pieces in \(z\)
        integer(ip) :: kq = 0_ip  !! The order of spline pieces in \(q\)
        integer(ip) :: kr = 0_ip  !! The order of spline pieces in \(r\)
        integer(ip) :: ks = 0_ip  !! The order of spline pieces in \(s\)
        real(wp),dimension(:,:,:,:,:,:),allocatable :: bcoef  !! array of coefficients of the b-spline interpolant
        real(wp),dimension(:),allocatable :: tx  !! The knots in the \(x\) direction for the spline interpolant
        real(wp),dimension(:),allocatable :: ty  !! The knots in the \(y\) direction for the spline interpolant
        real(wp),dimension(:),allocatable :: tz  !! The knots in the \(z\) direction for the spline interpolant
        real(wp),dimension(:),allocatable :: tq  !! The knots in the \(q\) direction for the spline interpolant
        real(wp),dimension(:),allocatable :: tr  !! The knots in the \(r\) direction for the spline interpolant
        real(wp),dimension(:),allocatable :: ts  !! The knots in the \(s\) direction for the spline interpolant
        integer(ip) :: inbvy = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: inbvz = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: inbvq = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: inbvr = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: inbvs = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: iloy  = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: iloz  = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: iloq  = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: ilor  = 1_ip  !! internal variable used for efficient processing
        integer(ip) :: ilos  = 1_ip  !! internal variable used for efficient processing
        real(wp),dimension(:,:,:,:,:),allocatable :: work_val_1  !! [[db6val]] work array of dimension `ky,kz,kq,kr,ks`
        real(wp),dimension(:,:,:,:),allocatable   :: work_val_2  !! [[db6val]] work array of dimension `kz,kq,kr,ks`
        real(wp),dimension(:,:,:),allocatable     :: work_val_3  !! [[db6val]] work array of dimension `kq,kr,ks`
        real(wp),dimension(:,:),allocatable       :: work_val_4  !! [[db6val]] work array of dimension `kr,ks`
        real(wp),dimension(:),allocatable         :: work_val_5  !! [[db6val]] work array of dimension `ks`
        real(wp),dimension(:),allocatable         :: work_val_6  !! [[db6val]] work array of dimension `3_ip*max(kx,ky,kz,kq,kr,ks)`
        contains
        private
        generic,public :: initialize => initialize_6d_auto_knots,initialize_6d_specify_knots
        procedure :: initialize_6d_auto_knots
        procedure :: initialize_6d_specify_knots
        procedure,public :: evaluate => evaluate_6d
        procedure,public :: destroy => destroy_6d
        procedure,public :: size_of => size_6d
        final :: finalize_6d
    end type bspline_6d