fitp_curvpd Function

private function fitp_curvpd(t, n, x, y, p, yp, sigma)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
real, intent(in) :: t
integer, intent(in) :: n
real, intent(in), dimension(:) :: x
real, intent(in), dimension(:) :: y
real, intent(in) :: p
real, intent(in), dimension(:) :: yp
real, intent(in) :: sigma

Return Value real


Contents

Source Code


Source Code

  real function fitp_curvpd (t,n,x,y,p,yp,sigma)
    real, intent(in) :: t, p, sigma
    integer, intent(in) :: n
    real, dimension(:), intent(in) :: x, y, yp
!
!                                 coded by alan kaylor cline
!                           from fitpack -- january 26, 1987
!                        a curve and surface fitting package
!                      a product of pleasant valley software
!                  8603 altus cove, austin, texas 78759, usa
!
! this function is the derivative of curvp2
! interpolates a curve at a given point using
! a periodic spline under tension. the subroutine curvp1
! should be called earlier to determine certain necessary
! parameters.
!
! on input--
!
!   t contains a real value to be mapped onto the interpo-
!   lating curve.
!
!   n contains the number of points which were specified to
!   determine the curve.
!
!   x and y are arrays containing the abscissae and
!   ordinates, respectively, of the specified points.
!
!   p contains the period.
!
!   yp is an array of second derivative values of the curve
!   at the nodes.
!
! and
!
!   sigma contains the tension factor (its sign is ignored).
!
! the parameters n, x, y, p, yp, and sigma should be input
! unaltered from the output of curvp1.
!
! on output--
!
!   curvpd contains the interpolated derivative
!
! none of the input parameters are altered.
!
! this function references package modules intrvp and
! snhcsh.
!
!-----------------------------------------------------------
    integer :: i, im1
    real :: ss, sigdel, sum, c2, c1, dummy
    real :: tp, sigmap, dels, del2, del1
!
! determine interval
!
    im1 = fitp_intrvp (t,x,n,p,tp)
    i = im1+1
!
! denormalize tension factor
!
    sigmap = abs(sigma)*real(n)/p
!
! set up and perform interpolation
!
    del1 = tp-x(im1)
    if (im1 .eq. n) go to 1
    del2 = x(i)-tp
    dels = x(i)-x(im1)
    go to 2
1   i = 1
    del2 = x(1)+p-tp
    dels = p-(x(n)-x(1))
2   sum = (y(i)-y(im1))/dels
    if (is_zero(sigmap)) then
      fitp_curvpd = sum+(yp(i)*(2.*del1*del1-del2*(del1+dels))- &
           yp(im1)*(2.*del2*del2-del1*(del2+dels))) &
           /(6.*dels)
    else
      sigdel = sigmap*dels
      call fitp_snhcsh (ss,dummy,sigdel,-1)
      call fitp_snhcsh (dummy,c1,sigmap*del1,1)
      call fitp_snhcsh (dummy,c2,sigmap*del2,1)
      fitp_curvpd = sum+(yp(i)*(c1-ss)-yp(im1)*(c2-ss))/(sigdel*sigmap*(1.+ss))
    end if
  end function fitp_curvpd