db1ink_alt_2 Subroutine

private pure subroutine db1ink_alt_2(x, nx, fcn, kx, ibcl, ibcr, fbcl, fbcr, tleft, tright, tx, bcoef, iflag)

Alternate version of db1ink_alt, where the first and last 3 knots are specified by the user.

History

  • Jacob Williams, 9/4/2018 : created this routine.

See also

  • dbint4 -- the main routine that is called here.

Note

Currently, this only works for 3rd order (k=4).

Arguments

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

vector of abscissae of length nx, distinct and in increasing order

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

number of data points,

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

vector of ordinates of length nx

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

spline order (Currently, this must be 4)

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

selection parameter for left boundary condition:

  • ibcl = 1 constrain the first derivative at x(1) to fbcl
  • ibcl = 2 constrain the second derivative at x(1) to fbcl
integer(kind=ip), intent(in) :: ibcr

selection parameter for right boundary condition:

  • ibcr = 1 constrain first derivative at x(nx) to fbcr
  • ibcr = 2 constrain second derivative at x(nx) to fbcr
real(kind=wp), intent(in) :: fbcl

left boundary values governed by ibcl

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

right boundary values governed by ibcr

real(kind=wp), intent(in), dimension(3) :: tleft

t(1:3) in increasing order supplied by the user.

real(kind=wp), intent(in), dimension(3) :: tright

t(nx+4:nx+6) in increasing order supplied by the user.

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

knot array of length nx+6

real(kind=wp), intent(out), dimension(:) :: bcoef

b spline coefficient array of length nx+2

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

status flag:

  • 0: no errors
  • 806: dbint4 can only be used when k=4

Calls

proc~~db1ink_alt_2~~CallsGraph proc~db1ink_alt_2 bspline_sub_module::db1ink_alt_2 proc~check_inputs bspline_sub_module::check_inputs proc~db1ink_alt_2->proc~check_inputs proc~dbint4 bspline_sub_module::dbint4 proc~db1ink_alt_2->proc~dbint4 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~dbspvn bspline_sub_module::dbspvn proc~dbspvd->proc~dbspvn

Called by

proc~~db1ink_alt_2~~CalledByGraph proc~db1ink_alt_2 bspline_sub_module::db1ink_alt_2 interface~db1ink bspline_sub_module::db1ink interface~db1ink->proc~db1ink_alt_2 proc~initialize_1d_auto_knots bspline_oo_module::bspline_1d%initialize_1d_auto_knots proc~initialize_1d_auto_knots->interface~db1ink proc~initialize_1d_specify_knots bspline_oo_module::bspline_1d%initialize_1d_specify_knots proc~initialize_1d_specify_knots->interface~db1ink proc~bspline_1d_constructor_auto_knots bspline_oo_module::bspline_1d_constructor_auto_knots proc~bspline_1d_constructor_auto_knots->proc~initialize_1d_auto_knots proc~bspline_1d_constructor_specify_knots bspline_oo_module::bspline_1d_constructor_specify_knots proc~bspline_1d_constructor_specify_knots->proc~initialize_1d_specify_knots interface~bspline_1d bspline_oo_module::bspline_1d interface~bspline_1d->proc~bspline_1d_constructor_auto_knots interface~bspline_1d->proc~bspline_1d_constructor_specify_knots

Source Code

    pure subroutine db1ink_alt_2(x,nx,fcn,kx,ibcl,ibcr,fbcl,fbcr,tleft,tright,tx,bcoef,iflag)

    implicit none

    real(wp),dimension(:),intent(in)   :: x       !! \(x\) vector of abscissae of length `nx`, distinct
                                                  !! and in increasing order
    integer(ip),intent(in)             :: nx      !! number of data points, \( n_x \ge 2 \)
    real(wp),dimension(:),intent(in)   :: fcn     !! \(y\) vector of ordinates of length `nx`
    integer(ip),intent(in)             :: kx      !! spline order (Currently, this must be `4`)
    integer(ip),intent(in)             :: ibcl    !! selection parameter for left boundary condition:
                                                  !!
                                                  !! * `ibcl = 1` constrain the first derivative at `x(1)` to `fbcl`
                                                  !! * `ibcl = 2` constrain the second derivative at `x(1)` to `fbcl`
    integer(ip),intent(in)             :: ibcr    !! selection parameter for right boundary condition:
                                                  !!
                                                  !! * `ibcr = 1` constrain first derivative at `x(nx)` to `fbcr`
                                                  !! * `ibcr = 2` constrain second derivative at `x(nx)` to `fbcr`
    real(wp),intent(in)                :: fbcl    !! left boundary values governed by `ibcl`
    real(wp),intent(in)                :: fbcr    !! right boundary values governed by `ibcr`
    real(wp),dimension(3),intent(in)   :: tleft   !! `t(1:3)` in increasing order supplied by the user.
    real(wp),dimension(3),intent(in)   :: tright  !! `t(nx+4:nx+6)` in increasing order supplied by the user.
    real(wp),dimension(:),intent(out)  :: tx      !! knot array of length `nx+6`
    real(wp),dimension(:),intent(out)  :: bcoef   !! b spline coefficient array of length `nx+2`
    integer(ip),intent(out)            :: iflag   !! status flag:
                                                  !!
                                                  !! * 0: no errors
                                                  !! * 806: [[dbint4]] can only be used when k=4

    real(wp),dimension(:,:),allocatable :: w         !! work array of dimension `5,nx+2`
    integer(ip)                         :: n         !! number of coefficients (`n=nx+2`)
    integer(ip)                         :: k         !! order of spline (`k=4`)
    logical                             :: status_ok !! status flag for error checking

    integer(ip),parameter :: kntopt = 3 !! use `tleft` and `tright` in [[dbint4]]

    if (kx /= 4_ip) then
        iflag = 806_ip
    else

        call check_inputs(  1_ip,& ! so it will check size of t
                            iflag,&
                            nx=nx,&
                            kx=kx,&
                            x=x,&
                            f1=fcn,&
                            bcoef1=bcoef,&
                            tx=tx,&
                            status_ok=status_ok,&
                            alt=.true.)

        if (status_ok) then
            allocate(w(5,nx+2))
            call dbint4(x,fcn,nx,ibcl,ibcr,fbcl,fbcr,kntopt,tleft,tright,tx,bcoef,n,k,w,iflag)
            deallocate(w)
        end if

    end if

    end subroutine db1ink_alt_2