Determines the parameters of a function that interpolates the two-dimensional gridded data The interpolating function and its derivatives may subsequently be evaluated by the function db2val.
The interpolating function is a piecewise polynomial function represented as a tensor product of one-dimensional b-splines. the form of this function is
where the functions and are one-dimensional b-spline basis functions. the coefficients are chosen so that
Note that for each fixed value of y, is a piecewise polynomial function of x alone, and for each fixed value of x, is a piecewise polynomial function of y alone. in one dimension a piecewise polynomial may be created by partitioning a given interval into subintervals and defining a distinct polynomial piece on each one. the points where adjacent subintervals meet are called knots. each of the functions and above is a piecewise polynomial.
Users of db2ink choose the order (degree+1) of the polynomial pieces used to define the piecewise polynomial in each of the x and y directions (kx and ky). users also may define their own knot sequence in x and y separately (tx and ty). if iflag=0, however, db2ink will choose sequences of knots that result in a piecewise polynomial interpolant with kx-2 continuous partial derivatives in x and ky-2 continuous partial derivatives in y. (kx knots are taken near each endpoint in the x direction, not-a-knot end conditions are used, and the remaining knots are placed at data points if kx is even or at midpoints between data points if kx is odd. the y direction is treated similarly.)
After a call to db2ink, all information necessary to define the interpolating function are contained in the parameters nx, ny, kx, ky, tx, ty, and bcoef. These quantities should not be altered until after the last call of the evaluation routine db2val.
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