Piecewise Cubic Hermite Function and Derivative evaluator
Evaluate a piecewise cubic Hermite function and its first derivative at an array of points. May be used by itself for Hermite interpolation, or as an evaluator for dpchim or dpchic.
Evaluates the cubic Hermite function defined by N, X, F, D, to- gether with its first derivative, at the points XE(J), J=1(1)NE.
If only function values are required, use dpchfe, instead.
To provide compatibility with dpchim and dpchic, includes an increment between successive values of the F- and D-arrays.
Most of the coding between the call to dchfdv and the end of the IR-loop could be eliminated if it were permissible to assume that XE is ordered relative to X.
dchfdv does not assume that X1 is less than X2. thus, it would be possible to write a version of DPCHFD that assumes a strictly decreasing X-array by simply running the IR-loop backwards (and reversing the order of appropriate tests).
The present code has a minor bug, which I have decided is not worth the effort that would be required to fix it. If XE contains points in [X(N-1),X(N)], followed by points < X(N-1), followed by points >X(N), the extrapolation points will be counted (at least) twice in the total returned in IERR.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
number of data points. (Error return if N<2). |
||
real(kind=wp), | intent(in) | :: | x(*) |
array of independent variable values. The
elements of X must be strictly increasing:
|
||
real(kind=wp), | intent(in) | :: | f(incfd,*) |
array of function values. F(1+(I-1)*INCFD) is the value corresponding to X(I). |
||
real(kind=wp), | intent(in) | :: | d(incfd,*) |
array of derivative values. D(1+(I-1)*INCFD) is the value corresponding to X(I). |
||
integer, | intent(in) | :: | incfd |
increment between successive values in F and D. (Error return if INCFD<1). |
||
logical, | intent(inout) | :: | skip |
logical variable which should be set to .TRUE. if the user wishes to skip checks for validity of preceding parameters, or to .FALSE. otherwise. This will save time in case these checks have already been performed (say, in dpchim or dpchic). SKIP will be set to .TRUE. on normal return. |
||
integer, | intent(in) | :: | ne |
number of evaluation points. (Error return if NE<1). |
||
real(kind=wp), | intent(in) | :: | xe(*) |
array of points at which the functions are to be evaluated. NOTES:
|
||
real(kind=wp), | intent(out) | :: | fe(*) |
array of values of the cubic Hermite function defined by N, X, F, D at the points XE. |
||
real(kind=wp), | intent(out) | :: | de(*) |
array of values of the first derivative of the same function at the points XE. |
||
integer, | intent(out) | :: | ierr |
error flag.
(Output arrays have not been changed in any of these cases.) NOTE: The above errors are checked in the order listed, and following arguments have NOT been validated.
|
subroutine dpchfd (n, x, f, d, incfd, skip, ne, xe, fe, de, ierr) integer,intent(in) :: n !! number of data points. (Error return if N<2). integer,intent(in) :: incfd !! increment between successive values in F and D. !! (Error return if INCFD<1). integer,intent(in) :: ne !! number of evaluation points. (Error return if NE<1). integer,intent(out) :: ierr !! error flag. !! !! * Normal return: !! * IERR = 0 (no errors). !! * Warning error: !! * IERR>0 means that extrapolation was performed at !! IERR points. !! * "Recoverable" errors: !! * IERR = -1 if N<2 . !! * IERR = -2 if INCFD<1 . !! * IERR = -3 if the X-array is not strictly increasing. !! * IERR = -4 if NE<1 . !! !! (Output arrays have not been changed in any of these cases.) !! !! NOTE: The above errors are checked in the order listed, !! and following arguments have **NOT** been validated. !! !! * IERR = -5 if an error has occurred in the lower-level !! routine DCHFDV. NB: this should never happen. !! Notify the author **IMMEDIATELY** if it does. real(wp),intent(in) :: x(*) !! array of independent variable values. The !! elements of X must be strictly increasing: !! `X(I-1) < X(I), I = 2(1)N`. !! (Error return if not.) real(wp),intent(in) :: f(incfd,*) !! array of function values. F(1+(I-1)*INCFD) is !! the value corresponding to X(I). real(wp),intent(in) :: d(incfd,*) !! array of derivative values. D(1+(I-1)*INCFD) !! is the value corresponding to X(I). real(wp),intent(in) :: xe(*) !! array of points at which the functions are to !! be evaluated. !! NOTES: !! !! 1. The evaluation will be most efficient if the elements !! of XE are increasing relative to X; !! that is, XE(J) >= X(I) !! implies XE(K) >= X(I), all K>=J . !! 2. If any of the XE are outside the interval [X(1),X(N)], !! values are extrapolated from the nearest extreme cubic, !! and a warning error is returned. real(wp),intent(out) :: fe(*) !! array of values of the cubic Hermite !! function defined by N, X, F, D at the points XE. real(wp),intent(out) :: de(*) !! array of values of the first derivative of !! the same function at the points XE. logical,intent(inout) :: skip !! logical variable which should be set to !! .TRUE. if the user wishes to skip checks for validity of !! preceding parameters, or to .FALSE. otherwise. !! This will save time in case these checks have already !! been performed (say, in [[dpchim]] or [[dpchic]]). !! SKIP will be set to .TRUE. on normal return. integer :: i, ierc, ir, j, jfirst, next(2), nj logical :: error, located ! validity-check arguments. if (.not. skip) then if ( n<2 ) then ierr = -1 call xermsg ('PCHIP', 'dpchfd', 'number of data points less than two', ierr, 1) return end if if ( incfd<1 ) then ierr = -2 call xermsg ('PCHIP', 'dpchfd', 'increment less than one', ierr, 1) return end if do i = 2, n if ( x(i)<=x(i-1) ) then ! x-array not strictly increasing. ierr = -3 call xermsg ('PCHIP', 'dpchfd', 'x-array not strictly increasing', ierr, 1) return end if end do end if ! function definition is ok, go on. if ( ne<1 ) then ierr = -4 call xermsg ('PCHIP', 'dpchfd', 'number of evaluation points less than one', ierr, 1) return end if ierr = 0 skip = .true. ! loop over intervals. ( interval index is il = ir-1 ) ! ( interval is x(il)<=x<x(ir) ) jfirst = 1 ir = 2 do ! skip out of loop if have processed all evaluation points. if (jfirst > ne) return ! locate all points in interval. located = .false. do j = jfirst, ne if (xe(j) >= x(ir)) then located = .true. exit end if end do if (located) then ! have located first point beyond interval. if (ir == n) j = ne + 1 else j = ne + 1 end if nj = j - jfirst ! skip evaluation if no points in interval. if (nj /= 0) then ! evaluate cubic at xe(i), i = jfirst (1) j-1 . ! ---------------------------------------------------------------- call dchfdv (x(ir-1),x(ir), f(1,ir-1),f(1,ir), d(1,ir-1), d(1,ir), nj, & xe(jfirst), fe(jfirst), de(jfirst), next, ierc) ! ---------------------------------------------------------------- if (ierc < 0) exit if (next(2) /= 0) then ! if (next(2) > 0) then ! in the current set of xe-points, there are next(2) to the ! right of x(ir). if (ir < n) exit ! if (ir == n) then ! these are actually extrapolation points. ierr = ierr + next(2) end if if (next(1) /= 0) then ! if (next(1) > 0) then ! in the current set of xe-points, there are next(1) to the ! left of x(ir-1). if (ir <= 2) then ! if (ir == 2) then ! these are actually extrapolation points. ierr = ierr + next(1) else ! xe is not ordered relative to x, so must adjust ! evaluation interval. ! first, locate first point to left of x(ir-1). error = .true. do i = jfirst, j-1 if (xe(i) < x(ir-1)) then error = .false. exit end if end do if (error) exit ! note-- cannot drop through here unless there is an error in dchfdv. ! reset j. (this will be the new jfirst.) j = i ! now find out how far to back up in the x-array. do i = 1, ir-1 if (xe(j) < x(i)) exit end do ! at this point, either ! xe(j) < x(1) or x(i-1) <= xe(j) < x(i) . ! reset ir, recognizing that it will be incremented before ! cycling. ir = max(1, i-1) end if end if jfirst = j ! end of ir-loop. end if ir = ir + 1 if (ir > n) return end do ! error return from dchfdv. ! *** this case should never occur *** ierr = -5 call xermsg ('PCHIP', 'dpchfd', 'error return from dchfdv -- fatal', ierr, 2) end subroutine dpchfd