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(ny) | :: | y | Array of y abcissae. Must be strictly increasing. |
|
integer, | intent(in) | :: | ny | Number of y abcissae |
||
real(kind=wp), | intent(in), | dimension(nx,ny) | :: | fcn | Array of function values to interpolate. fcn(i,j) should contain the function value at the point (x(i),y(j)) |
|
integer, | intent(in) | :: | kx | The order of spline pieces in x (>= 2, < nx). (order = polynomial degree + 1) |
||
integer, | intent(in) | :: | ky | The order of spline pieces in y (>= 2, < ny). (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 db2ink. If iflag=1 these are specified by the user. Must be non-decreasing. |
|
real(kind=wp), | intent(inout), | dimension(ny+ky) | :: | ty | The knots in the y direction for the spline interpolant. If iflag=0 these are chosen by db2ink. If iflag=1 these are specified by the user. Must be non-decreasing. |
|
real(kind=wp), | intent(out), | dimension(nx,ny) | :: | bcoef | Array of coefficients of the b-spline interpolant. |
|
integer, | intent(inout) | :: | iflag | on input: 0 = knot sequence chosen by db2ink. 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. 7 = ny out of range. 8 = ky out of range. 9 = y not strictly increasing. 10 = ty not non-decreasing. |
subroutine db2ink(x,nx,y,ny,fcn,kx,ky,tx,ty,bcoef,iflag)
implicit none
integer,intent(in) :: nx !! Number of x abcissae
integer,intent(in) :: ny !! Number of y abcissae
integer,intent(in) :: kx !! The order of spline pieces in x (>= 2, < nx). (order = polynomial degree + 1)
integer,intent(in) :: ky !! The order of spline pieces in y (>= 2, < ny). (order = polynomial degree + 1)
real(wp),dimension(nx),intent(in) :: x !! Array of x abcissae. Must be strictly increasing.
real(wp),dimension(ny),intent(in) :: y !! Array of y abcissae. Must be strictly increasing.
real(wp),dimension(nx,ny),intent(in) :: fcn !! Array of function values to interpolate. fcn(i,j) should
!! contain the function value at the point (x(i),y(j))
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 [[db2ink]].
!! If iflag=1 these are specified by the user.
!! Must be non-decreasing.
real(wp),dimension(ny+ky),intent(inout) :: ty !! The knots in the y direction for the spline interpolant.
!! If iflag=0 these are chosen by [[db2ink]].
!! If iflag=1 these are specified by the user.
!! Must be non-decreasing.
real(wp),dimension(nx,ny),intent(out) :: bcoef !! Array of coefficients of the b-spline interpolant.
integer,intent(inout) :: iflag !! **on input:** 0 = knot sequence chosen by [[db2ink]].
!! 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.
!! 7 = ny out of range.
!! 8 = ky out of range.
!! 9 = y not strictly increasing.
!! 10 = ty not non-decreasing.
real(wp),dimension(nx*ny) :: temp
real(wp),dimension(max(2*kx*(nx+1),2*ky*(ny+1))) :: work
logical :: status_ok
!check validity of inputs
call check_inputs('db2ink',&
iflag,&
nx=nx,ny=ny,&
kx=kx,ky=ky,&
x=x,y=y,&
tx=tx,ty=ty,&
status_ok=status_ok)
if (status_ok) then
!choose knots
if (iflag == 0) then
call dbknot(x,nx,kx,tx)
call dbknot(y,ny,ky,ty)
endif
!construct b-spline coefficients
call dbtpcf(x,nx,fcn, nx,ny,tx,kx,temp, work,iflag)
if (iflag==0) call dbtpcf(y,ny,temp,ny,nx,ty,ky,bcoef,work,iflag)
if (iflag==0) iflag = 1
endif
endsubroutine db2ink