UltrasphericalGradientEvalSum
Evaluate finite sum of gradient of Ultraspherical polynomials.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalGradientEvalSum(n, lambda, x, &
& coeff) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda 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 UltrasphericalGradientEvalSum
END INTERFACE
This example shows the usage of UltrasphericalGradientEvalSum
method.
This routine evaluates the gradient of 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( n=2_I4B, x=x, &
& lambda=lambda ) &
& + 2.0_DFP * UltrasphericalGradientEval( &
& 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 UltrasphericalGradientEvalSum(n, lambda, x, coeff) &
& RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda 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 UltrasphericalGradientEvalSum
END INTERFACE
This example shows the usage of UltrasphericalGradientEvalSum
method.
This routine evaluates the gradient offinite 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff )
exact = UltrasphericalGradientEval( n=n, x=x, &
& lambda=lambda )
call ok( ALL(SOFTEQ(ans, exact, tol )))
end program main
Interface 3
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalGradientEvalSum(n, lambda, x, &
& coeff, k) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda of Ultraspherical polynomial
REAL(DFP), INTENT(IN) :: x
!! point
REAL(DFP), INTENT(IN) :: coeff(0:n)
!! Coefficient of finite sum, size = n+1
INTEGER(I4B), INTENT(IN) :: k
!! order of derivative
REAL(DFP) :: ans
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalGradientEvalSum
END INTERFACE
This example shows the usage of UltrasphericalGradientEvalSum
method.
This routine evaluates the kth gradient of 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff, k=1 )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff, k=1 )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff, k=1 )
exact = UltrasphericalGradientEval( 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 = UltrasphericalGradientEvalSum( n=n, x=x, &
& lambda=lambda, coeff=coeff, k = 1 )
exact = UltrasphericalGradientEval( n=2_I4B, x=x, &
& lambda=lambda ) &
& + 2.0_DFP * UltrasphericalGradientEval( &
& n=3_I4B, x=x, &
& lambda=lambda )
call ok( SOFTEQ(ans, exact, tol ))
end program main
Interface 4
INTERFACE
MODULE PURE FUNCTION UltrasphericalGradientEvalSum(n, lambda, x, &
& coeff, k) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda of Ultraspherical polynomial
REAL(DFP), INTENT(IN) :: x(:)
!! point
REAL(DFP), INTENT(IN) :: coeff(0:n)
!! Coefficient of finite sum, size = n+1
INTEGER(I4B), INTENT(IN) :: k
!! kth order derivative
REAL(DFP) :: ans(SIZE(x))
!! Evaluate Ultraspherical polynomial of order n at point x
END FUNCTION UltrasphericalGradientEvalSum
END INTERFACE