Initialize a bspline_2d type (with user-specified knots). This is a wrapper for db2ink.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(bspline_2d), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in), | dimension(:) | :: | x |
|
|
real(kind=wp), | intent(in), | dimension(:) | :: | y |
|
|
real(kind=wp), | intent(in), | dimension(:,:) | :: | fcn |
|
|
integer(kind=ip), | intent(in) | :: | kx |
The order of spline pieces in ( ) (order = polynomial degree + 1) |
||
integer(kind=ip), | intent(in) | :: | ky |
The order of spline pieces in ( ) (order = polynomial degree + 1) |
||
real(kind=wp), | intent(in), | dimension(:) | :: | tx |
The |
|
real(kind=wp), | intent(in), | dimension(:) | :: | ty |
The |
|
integer(kind=ip), | intent(out) | :: | iflag |
status flag (see db2ink) |
||
logical, | intent(in), | optional | :: | extrap |
if true, then extrapolation is allowed (default is false) |
pure subroutine initialize_2d_specify_knots(me,x,y,fcn,kx,ky,tx,ty,iflag,extrap) implicit none class(bspline_2d),intent(inout) :: me real(wp),dimension(:),intent(in) :: x !! `(nx)` array of \(x\) abcissae. Must be strictly increasing. real(wp),dimension(:),intent(in) :: y !! `(ny)` array of \(y\) abcissae. Must be strictly increasing. real(wp),dimension(:,:),intent(in) :: fcn !! `(nx,ny)` matrix of function values to interpolate. !! `fcn(i,j)` should contain the function value at the !! point (`x(i)`,`y(j)`) integer(ip),intent(in) :: kx !! The order of spline pieces in \(x\) !! ( \( 2 \le k_x < n_x \) ) !! (order = polynomial degree + 1) integer(ip),intent(in) :: ky !! The order of spline pieces in \(y\) !! ( \( 2 \le k_y < n_y \) ) !! (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. real(wp),dimension(:),intent(in) :: ty !! The `(ny+ky)` knots in the \(y\) direction !! for the spline interpolant. !! Must be non-decreasing. integer(ip),intent(out) :: iflag !! status flag (see [[db2ink]]) logical,intent(in),optional :: extrap !! if true, then extrapolation is allowed !! (default is false) integer(ip) :: nx,ny call me%destroy() nx = size(x,kind=ip) ny = size(y,kind=ip) call check_knot_vectors_sizes(nx=nx,kx=kx,tx=tx,& ny=ny,ky=ky,ty=ty,& iflag=iflag) if (iflag == 0_ip) then me%nx = nx me%ny = ny me%kx = kx me%ky = ky allocate(me%tx(nx+kx)) allocate(me%ty(ny+ky)) allocate(me%bcoef(nx,ny)) allocate(me%work_val_1(ky)) allocate(me%work_val_2(3_ip*max(kx,ky))) me%tx = tx me%ty = ty call db2ink(x,nx,y,ny,fcn,kx,ky,1_ip,me%tx,me%ty,me%bcoef,iflag) call me%set_extrap_flag(extrap) end if me%initialized = iflag==0_ip me%iflag = iflag end subroutine initialize_2d_specify_knots