UltrasphericalEval
Evaluate Ultraspherical polynomials of order n at single or several points.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEval(n, lambda, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda should be greater than -0.5
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalEval
END INTERFACE
Evaluate Ultraspherical polynomial of order n at single points.
- N, the order of polynomial to compute.
- lambda is the polynomial parameter.
- x: the point at which the polynomials are to be evaluated.
This example shows the usage of UltrasphericalEval
method.
This routine evaluates Ultraspherical polynomial of order n, at single point
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ) :: ans, x, exact
real( dfp ), parameter :: tol=1.0E-10, lambda=0.5_DFP
n = 3
x = -1.0_DFP; call callme
exact = 0.5_DFP*(5.0 * x**3 - 3.0*x)
call ok( SOFTEQ(ans, exact, tol ))
x = 0.0_DFP; call callme
exact = 0.5_DFP*(5.0 * x**3 - 3.0*x)
call ok( SOFTEQ(ans, exact, tol ))
x = +1.0_DFP; call callme
exact = 0.5_DFP*(5.0 * x**3 - 3.0*x)
call ok( SOFTEQ(ans, exact, tol ))
contains
subroutine callme
ans= UltrasphericalEval( n=n, x=x, lambda=lambda )
end subroutine callme
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEval(n, lambda, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda should be greater than -0.5
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(SIZE(x))
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalEval
END INTERFACE
- N is order of polynomial to compute.
- lambda is the polynomial parameter.
- x: the point at which the polynomials are to be evaluated.
- ans, the values of the Ultraspherical polynomials at the several points.
This example shows the usage of UltrasphericalEval
method.
This routine evaluates Ultraspherical polynomial of order n, at several points.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans(:), x(:), exact(:)
real( dfp ), parameter :: tol=1.0E-10, lambda=0.5_DFP
type(string) :: astr
n = 3
x = [-1.0, -0.5, 0.0, 0.5, 1.0]; call callme
exact = 0.5_DFP*(5.0 * x**3 - 3.0*x)
call ok( ALL(SOFTEQ(ans, exact, tol )))
contains
subroutine callme
ans= UltrasphericalEval( n=n, x=x, lambda=lambda )
end subroutine callme
end program main