EZspline_interp


interface EZspline_interp
!
! Interpolation at grid point(s) p1, [p2, [p3]]. Result is returned in
! f. Interpolation can be sought at a single point (p1, [p2, [p3]] are
! scalars), on an unordered list of points (p1, [p2, [p3]] have dimension
! k), or on a structured grid (p1, [p2, [p3]] have dimension k1, [k2, [k3]]
! respectively).
!
subroutine EZspline_interp3_r8(spline_o, p1, p2, p3, f, ier)
! single point evaluation
use EZspline_obj
type(EZspline3_r8) spline_o
real(ezspline_r8) :: p1, p2, p3
real(ezspline_r8) f
integer, intent(out) :: ier
end subroutine EZspline_interp3_r8

subroutine EZspline_interp3_array_r8(spline_o, k1, k2, k3, p1, p2, p3, f, ier)
use EZspline_obj
type(EZspline3_r8) spline_o
integer :: k1, k2, k3
real(ezspline_r8), intent(in) :: p1(k1), p2(k2), p3(k3)
real(ezspline_r8), intent(out):: f(k1,k2,*)
integer, intent(out) :: ier
end subroutine EZspline_interp3_array_r8

subroutine EZspline_interp3_cloud_r8(spline_o, k, p1, p2, p3, f, ier)
use EZspline_obj
type(EZspline3_r8) spline_o
integer, intent(in) :: k
real(ezspline_r8), intent(in) :: p1(k), p2(k), p3(k)
real(ezspline_r8), intent(out):: f(k)
integer, intent(out) :: ier
end subroutine EZspline_interp3_cloud_r8

subroutine EZspline_interp2_r8(spline_o, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r8) spline_o
real(ezspline_r8) :: p1, p2
real(ezspline_r8) f
integer, intent(out) :: ier
end subroutine EZspline_interp2_r8

subroutine EZspline_interp2_array_r8(spline_o, k1, k2, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r8) spline_o
integer :: k1, k2
real(ezspline_r8), intent(in) :: p1(k1), p2(k2)
real(ezspline_r8), intent(out):: f(k1,*)
integer, intent(out) :: ier
end subroutine EZspline_interp2_array_r8

subroutine EZspline_interp2_cloud_r8(spline_o, k, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r8) spline_o
integer, intent(in) :: k
real(ezspline_r8), intent(in) :: p1(k), p2(k)
real(ezspline_r8), intent(out):: f(k)
integer, intent(out) :: ier
end subroutine EZspline_interp2_cloud_r8

subroutine EZspline_interp1_r8(spline_o, p1, f, ier)
use EZspline_obj
type(EZspline1_r8) spline_o
real(ezspline_r8) :: p1
real(ezspline_r8) f
integer, intent(out) :: ier
end subroutine EZspline_interp1_r8

subroutine EZspline_interp1_array_r8(spline_o, k1, p1, f, ier)
use EZspline_obj
type(EZspline1_r8) spline_o
integer :: k1
real(ezspline_r8), intent(in) :: p1(k1)
real(ezspline_r8), intent(out):: f(k1)
integer, intent(out) :: ier
end subroutine EZspline_interp1_array_r8

subroutine EZspline_interp3_r4(spline_o, p1, p2, p3, f, ier)
! single point evaluation
use EZspline_obj
type(EZspline3_r4) spline_o
real(ezspline_r4) :: p1, p2, p3
real(ezspline_r4) f
integer, intent(out) :: ier
end subroutine EZspline_interp3_r4

subroutine EZspline_interp3_array_r4(spline_o, k1, k2, k3, p1, p2, p3, f, ier)
use EZspline_obj
type(EZspline3_r4) spline_o
integer :: k1, k2, k3
real(ezspline_r4), intent(in) :: p1(k1), p2(k2), p3(k3)
real(ezspline_r4), intent(out):: f(k1,k2,*)
integer, intent(out) :: ier
end subroutine EZspline_interp3_array_r4

subroutine EZspline_interp3_cloud_r4(spline_o, k, p1, p2, p3, f, ier)
use EZspline_obj
type(EZspline3_r4) spline_o
integer, intent(in) :: k
real(ezspline_r4), intent(in) :: p1(k), p2(k), p3(k)
real(ezspline_r4), intent(out):: f(k)
integer, intent(out) :: ier
end subroutine EZspline_interp3_cloud_r4

subroutine EZspline_interp2_r4(spline_o, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r4) spline_o
real(ezspline_r4) :: p1, p2
real(ezspline_r4) f
integer, intent(out) :: ier
end subroutine EZspline_interp2_r4

subroutine EZspline_interp2_array_r4(spline_o, k1, k2, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r4) spline_o
integer :: k1, k2
real(ezspline_r4), intent(in) :: p1(k1), p2(k2)
real(ezspline_r4), intent(out):: f(k1,*)
integer, intent(out) :: ier
end subroutine EZspline_interp2_array_r4

subroutine EZspline_interp2_cloud_r4(spline_o, k, p1, p2, f, ier)
use EZspline_obj
type(EZspline2_r4) spline_o
integer, intent(in) :: k
real(ezspline_r4), intent(in) :: p1(k), p2(k)
real(ezspline_r4), intent(out):: f(k)
integer, intent(out) :: ier
end subroutine EZspline_interp2_cloud_r4

subroutine EZspline_interp1_r4(spline_o, p1, f, ier)
use EZspline_obj
type(EZspline1_r4) spline_o
real(ezspline_r4) :: p1
real(ezspline_r4) f
integer, intent(out) :: ier
end subroutine EZspline_interp1_r4

subroutine EZspline_interp1_array_r4(spline_o, k1, p1, f, ier)
use EZspline_obj
type(EZspline1_r4) spline_o
integer :: k1
real(ezspline_r4), intent(in) :: p1(k1)
real(ezspline_r4), intent(out):: f(k1)
integer, intent(out) :: ier
end subroutine EZspline_interp1_array_r4

end interface


Send comments about this document to pletzer@pppl.gov. Tue Apr 24 14:01:03 2007