Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(nx) | :: | x | Array of x abcissae. Must be strictly increasing. |
|
integer, | intent(in) | :: | nx | Number of x abcissae |
||
real(kind=wp), | intent(in), | dimension(nx) | :: | fcn | Array of function values to interpolate. fcn(i) should contain the function value at the point x(i) |
|
integer, | intent(in) | :: | kx | The order of spline pieces in x (>= 2, < nx). (order = polynomial degree + 1) |
||
real(kind=wp), | intent(inout), | dimension(nx+kx) | :: | tx | The knots in the x direction for the spline interpolant. If iflag=0 these are chosen by db1ink. If iflag=1 these are specified by the user. Must be non-decreasing. |
|
real(kind=wp), | intent(out), | dimension(nx) | :: | bcoef | Array of coefficients of the b-spline interpolant. |
|
integer, | intent(inout) | :: | iflag | on input: 0 = knot sequence chosen by db1ink. 1 = knot sequence chosen by user. on output: 1 = successful execution. 2 = iflag out of range. 3 = nx out of range. 4 = kx out of range. 5 = x not strictly increasing. 6 = tx not non-decreasing. |
subroutine db1ink(x,nx,fcn,kx,tx,bcoef,iflag)
implicit none
integer,intent(in) :: nx !! Number of x abcissae
integer,intent(in) :: kx !! The order of spline pieces in x (>= 2, < nx). (order = polynomial degree + 1)
real(wp),dimension(nx),intent(in) :: x !! Array of x abcissae. Must be strictly increasing.
real(wp),dimension(nx),intent(in) :: fcn !! Array of function values to interpolate. fcn(i) should
!! contain the function value at the point x(i)
real(wp),dimension(nx+kx),intent(inout) :: tx !! The knots in the x direction for the spline interpolant.
!! If iflag=0 these are chosen by [[db1ink]].
!! If iflag=1 these are specified by the user.
!! Must be non-decreasing.
real(wp),dimension(nx),intent(out) :: bcoef !! Array of coefficients of the b-spline interpolant.
integer,intent(inout) :: iflag !! **on input:** 0 = knot sequence chosen by [[db1ink]].
!! 1 = knot sequence chosen by user.
!! **on output:** 1 = successful execution.
!! 2 = iflag out of range.
!! 3 = nx out of range.
!! 4 = kx out of range.
!! 5 = x not strictly increasing.
!! 6 = tx not non-decreasing.
real(wp),dimension(2*kx*(nx+1)) :: work
logical :: status_ok
!check validity of inputs
call check_inputs('db1ink',&
iflag,&
nx=nx,&
kx=kx,&
x=x,&
tx=tx,&
status_ok=status_ok)
if (status_ok) then
!choose knots
if (iflag == 0) then
call dbknot(x,nx,kx,tx)
endif
!construct b-spline coefficients
call dbtpcf(x,nx,fcn,nx,1,tx,kx,bcoef,work,iflag)
if (iflag==0) iflag = 1
endif
endsubroutine db1ink