Determines the parameters of a function that interpolates the one-dimensional gridded data The interpolating function and its derivatives may subsequently be evaluated by the function db1val.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(:) | :: | x |
|
|
integer(kind=ip), | intent(in) | :: | nx |
Number of abcissae |
||
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) | :: | iknot |
knot sequence flag:
|
||
real(kind=wp), | intent(inout), | dimension(:) | :: | tx |
The
Must be non-decreasing. |
|
real(kind=wp), | intent(out), | dimension(:) | :: | bcoef |
|
|
integer(kind=ip), | intent(out) | :: | iflag |
status flag:
|
pure subroutine db1ink_default(x,nx,fcn,kx,iknot,tx,bcoef,iflag) implicit none integer(ip),intent(in) :: nx !! Number of \(x\) abcissae integer(ip),intent(in) :: kx !! The order of spline pieces in \(x\) !! ( \( 2 \le k_x < n_x \) ) !! (order = polynomial degree + 1) real(wp),dimension(:),intent(in) :: x !! `(nx)` array of \(x\) abcissae. Must be strictly increasing. real(wp),dimension(:),intent(in) :: fcn !! `(nx)` array of function values to interpolate. `fcn(i)` should !! contain the function value at the point `x(i)` integer(ip),intent(in) :: iknot !! knot sequence flag: !! !! * 0 = knot sequence chosen by [[db1ink]]. !! * 1 = knot sequence chosen by user. real(wp),dimension(:),intent(inout) :: tx !! The `(nx+kx)` knots in the \(x\) direction !! for the spline interpolant: !! !! * If `iknot=0` these are chosen by [[db1ink]]. !! * If `iknot=1` these are specified by the user. !! !! Must be non-decreasing. real(wp),dimension(:),intent(out) :: bcoef !! `(nx)` array of coefficients of the b-spline interpolant. integer(ip),intent(out) :: iflag !! status flag: !! !! * 0 = successful execution. !! * 2 = `iknot` out of range. !! * 3 = `nx` out of range. !! * 4 = `kx` out of range. !! * 5 = `x` not strictly increasing. !! * 6 = `tx` not non-decreasing. !! * 700 = `size(x)` \( \ne \) `size(fcn,1)`. !! * 706 = `size(x)` \( \ne \) `nx`. !! * 712 = `size(tx)` \( \ne \) `nx+kx`. !! * 800 = `size(x)` \( \ne \) `size(bcoef,1)`. logical :: status_ok real(wp),dimension(:),allocatable :: work !! work array of dimension `2*kx*(nx+1)` !check validity of inputs call check_inputs( iknot,& iflag,& nx=nx,& kx=kx,& x=x,& f1=fcn,& bcoef1=bcoef,& tx=tx,& status_ok=status_ok) if (status_ok) then !choose knots if (iknot == 0_ip) then call dbknot(x,nx,kx,tx) end if allocate(work(2_ip*kx*(nx+1_ip))) !construct b-spline coefficients call dbtpcf(x,nx,fcn,nx,1_ip,tx,kx,bcoef,work,iflag) deallocate(work) end if end subroutine db1ink_default