UltrasphericalEvalSum
Evaluate the finite sum of Ultraspherical polynomials at point x.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEvalSum(n, lambda, x, coeff) &
& RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! alpha of Ultraspherical polynomial
REAL(DFP), INTENT(IN) :: x
!! point
REAL(DFP), INTENT(IN) :: coeff(0:n)
!! Coefficient of finite sum, size = n+1
REAL(DFP) :: ans
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalEvalSum
END INTERFACE
This example shows the usage of UltrasphericalEvalSum
method.
This routine evaluates sum of finite series of Ultraspherical polynomials of order upto n, at a single point.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ) :: ans, x, exact
real( dfp ), allocatable :: coeff(:)
real( dfp ), parameter :: tol=1.0E-10, lambda=0.5_DFP
n = 3; coeff = [1, 0, 0, 0]
x = -0.5_DFP
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=0_I4B, x=x, &
& lambda=lambda )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 1, 0, 0]
x = -0.5_DFP
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=1_I4B, x=x, &
& lambda=lambda )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 0, 1]
x = 0.25_DFP
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=3_I4B, x=x, &
& lambda=lambda )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 1, 2]
x = -0.5_DFP
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=2_I4B, x=x, &
& lambda=lambda ) &
& + 2.0_DFP * UltrasphericalEval( n=3_I4B, x=x, &
& lambda=lambda )
call ok( SOFTEQ(ans, exact, tol ))
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEvalSum(n, lambda, x, coeff) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! alpha of Ultraspherical polynomial
REAL(DFP), INTENT(IN) :: x(:)
!! point
REAL(DFP), INTENT(IN) :: coeff(0:n)
!! Coefficient of finite sum, size = n+1
REAL(DFP) :: ans(SIZE(x))
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalEvalSum
END INTERFACE
This example shows the usage of UltrasphericalEvalSum
method.
This routine evaluates the finite sum of the Ultraspherical polynomials of order upto n, at several points.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans(:), x(:), exact(:), &
& coeff(:)
real( dfp ), parameter :: tol=1.0E-10, lambda=0.5_DFP
type(string) :: astr
n = 3
x = [-0.5, 0.5, 0.25]
coeff = [1, 0, 0, 0]
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=0_I4B, x=x, &
& lambda=lambda )
call ok( ALL(SOFTEQ(ans, exact, tol )))
n = 3
x = [-0.5, 0.5, 0.25]
coeff = [0, 0, 0, 1]
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=3_I4B, x=x, &
& lambda=lambda )
call ok( ALL(SOFTEQ(ans, exact, tol )))
n = 100
x = [-0.5, 0.5, 0.25]
coeff = zeros(n+1, 0.0_DFP)
coeff( n+1 ) = 1.0_DFP
ans = UltrasphericalEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalEval( n=n, x=x, &
& lambda=lambda )
call ok( ALL(SOFTEQ(ans, exact, tol )))
end program main