bspline_1d_constructor_specify_knots Function

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)


Calls

proc~~bspline_1d_constructor_specify_knots~~CallsGraph proc~bspline_1d_constructor_specify_knots bspline_oo_module::bspline_1d_constructor_specify_knots proc~initialize_1d_specify_knots bspline_oo_module::bspline_1d%initialize_1d_specify_knots proc~bspline_1d_constructor_specify_knots->proc~initialize_1d_specify_knots interface~db1ink bspline_sub_module::db1ink proc~initialize_1d_specify_knots->interface~db1ink proc~check_knot_vectors_sizes bspline_oo_module::check_knot_vectors_sizes proc~initialize_1d_specify_knots->proc~check_knot_vectors_sizes proc~destroy_1d bspline_oo_module::bspline_1d%destroy_1d proc~initialize_1d_specify_knots->proc~destroy_1d proc~set_extrap_flag bspline_oo_module::bspline_class%set_extrap_flag proc~initialize_1d_specify_knots->proc~set_extrap_flag proc~db1ink_alt bspline_sub_module::db1ink_alt interface~db1ink->proc~db1ink_alt proc~db1ink_alt_2 bspline_sub_module::db1ink_alt_2 interface~db1ink->proc~db1ink_alt_2 proc~db1ink_default bspline_sub_module::db1ink_default interface~db1ink->proc~db1ink_default proc~destroy_base bspline_oo_module::bspline_class%destroy_base proc~destroy_1d->proc~destroy_base proc~check_inputs bspline_sub_module::check_inputs proc~db1ink_alt->proc~check_inputs proc~dbint4 bspline_sub_module::dbint4 proc~db1ink_alt->proc~dbint4 proc~db1ink_alt_2->proc~check_inputs proc~db1ink_alt_2->proc~dbint4 proc~db1ink_default->proc~check_inputs proc~dbknot bspline_sub_module::dbknot proc~db1ink_default->proc~dbknot proc~dbtpcf bspline_sub_module::dbtpcf proc~db1ink_default->proc~dbtpcf proc~dbnfac bspline_sub_module::dbnfac proc~dbint4->proc~dbnfac proc~dbnslv bspline_sub_module::dbnslv proc~dbint4->proc~dbnslv proc~dbspvd bspline_sub_module::dbspvd proc~dbint4->proc~dbspvd proc~dbintk bspline_sub_module::dbintk proc~dbtpcf->proc~dbintk proc~dbtpcf->proc~dbnslv proc~dbintk->proc~dbnfac proc~dbintk->proc~dbnslv proc~dbspvn bspline_sub_module::dbspvn proc~dbintk->proc~dbspvn proc~dbspvd->proc~dbspvn

Called by

proc~~bspline_1d_constructor_specify_knots~~CalledByGraph proc~bspline_1d_constructor_specify_knots bspline_oo_module::bspline_1d_constructor_specify_knots interface~bspline_1d bspline_oo_module::bspline_1d interface~bspline_1d->proc~bspline_1d_constructor_specify_knots

Source Code

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

    implicit none

    type(bspline_1d)                 :: me
    real(wp),dimension(:),intent(in) :: x     !! `(nx)` array of \(x\) abcissae. Must be strictly increasing.
    real(wp),dimension(:),intent(in) :: fcn   !! `(nx)` array of function values to interpolate. `fcn(i)` should
                                              !! contain the function value at the point `x(i)`
    integer(ip),intent(in)           :: kx    !! The order of spline pieces in \(x\)
                                              !! ( \( 2 \le k_x < n_x \) )
                                              !! (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.
    logical,intent(in),optional      :: extrap !! if true, then extrapolation is allowed
                                               !! (default is false)

    call initialize_1d_specify_knots(me,x,fcn,kx,tx,me%iflag,extrap)

    end function bspline_1d_constructor_specify_knots