bspline_2d_constructor_specify_knots Function

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)


Calls

proc~~bspline_2d_constructor_specify_knots~~CallsGraph proc~bspline_2d_constructor_specify_knots bspline_oo_module::bspline_2d_constructor_specify_knots proc~initialize_2d_specify_knots bspline_oo_module::bspline_2d%initialize_2d_specify_knots proc~bspline_2d_constructor_specify_knots->proc~initialize_2d_specify_knots proc~check_knot_vectors_sizes bspline_oo_module::check_knot_vectors_sizes proc~initialize_2d_specify_knots->proc~check_knot_vectors_sizes proc~db2ink bspline_sub_module::db2ink proc~initialize_2d_specify_knots->proc~db2ink proc~destroy_2d bspline_oo_module::bspline_2d%destroy_2d proc~initialize_2d_specify_knots->proc~destroy_2d proc~set_extrap_flag bspline_oo_module::bspline_class%set_extrap_flag proc~initialize_2d_specify_knots->proc~set_extrap_flag proc~check_inputs bspline_sub_module::check_inputs proc~db2ink->proc~check_inputs proc~dbknot bspline_sub_module::dbknot proc~db2ink->proc~dbknot proc~dbtpcf bspline_sub_module::dbtpcf proc~db2ink->proc~dbtpcf proc~destroy_base bspline_oo_module::bspline_class%destroy_base proc~destroy_2d->proc~destroy_base proc~dbintk bspline_sub_module::dbintk proc~dbtpcf->proc~dbintk proc~dbnslv bspline_sub_module::dbnslv proc~dbtpcf->proc~dbnslv proc~dbintk->proc~dbnslv proc~dbnfac bspline_sub_module::dbnfac proc~dbintk->proc~dbnfac proc~dbspvn bspline_sub_module::dbspvn proc~dbintk->proc~dbspvn

Called by

proc~~bspline_2d_constructor_specify_knots~~CalledByGraph proc~bspline_2d_constructor_specify_knots bspline_oo_module::bspline_2d_constructor_specify_knots interface~bspline_2d bspline_oo_module::bspline_2d interface~bspline_2d->proc~bspline_2d_constructor_specify_knots

Source Code

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

    implicit none

    type(bspline_2d)                   :: me
    real(wp),dimension(:),intent(in)   :: x     !! `(nx)` array of \(x\) abcissae. Must be strictly increasing.
    real(wp),dimension(:),intent(in)   :: y     !! `(ny)` array of \(y\) abcissae. Must be strictly increasing.
    real(wp),dimension(:,:),intent(in) :: 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(ip),intent(in)             :: kx    !! The order of spline pieces in \(x\)
                                                !! ( \( 2 \le k_x < n_x \) )
                                                !! (order = polynomial degree + 1)
    integer(ip),intent(in)             :: ky    !! The order of spline pieces in \(y\)
                                                !! ( \( 2 \le k_y < n_y \) )
                                                !! (order = polynomial degree + 1)
    real(wp),dimension(:),intent(in)   :: tx    !! The `(nx+kx)` knots in the \(x\) direction
                                                !! for the spline interpolant.
                                                !! Must be non-decreasing.
    real(wp),dimension(:),intent(in)   :: ty    !! The `(ny+ky)` knots in the \(y\) direction
                                                !! for the spline interpolant.
                                                !! Must be non-decreasing.
    logical,intent(in),optional      :: extrap  !! if true, then extrapolation is allowed
                                                !! (default is false)

    call initialize_2d_specify_knots(me,x,y,fcn,kx,ky,tx,ty,me%iflag,extrap)

    end function bspline_2d_constructor_specify_knots