Chebyshev1GradientEvalSum
Evaluate finite sum of gradient of Chebyshev1 polynomials.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1GradientEvalSum(n, x, coeff) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of 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 Chebyshev1 polynomial of order n at point x
END FUNCTION Chebyshev1GradientEvalSum
END INTERFACE
This example shows the usage of Chebyshev1GradientEvalSum
method.
This routine evaluates the gradient of sum of finite series of Chebyshev1 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
n = 3; coeff = [1, 0, 0, 0]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=0_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 1, 0, 0]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=1_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 0, 1]
x = 0.25_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=3_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 1, 2]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=2_I4B, x=x) &
& + 2.0_DFP * Chebyshev1GradientEval( &
& n=3_I4B, x=x)
call ok( SOFTEQ(ans, exact, tol ))
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1GradientEvalSum(n, x, coeff) &
& RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of 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 Chebyshev1 polynomial of order n at point x
END FUNCTION Chebyshev1GradientEvalSum
END INTERFACE
This example shows the usage of Chebyshev1GradientEvalSum
method.
This routine evaluates the gradient offinite sum of the Chebyshev1 polynomials of order upto n, at several point.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans(:), x(:), exact(:), &
& coeff(:)
real( dfp ), parameter :: tol=1.0E-10
type(string) :: astr
n = 3
x = [-0.5, 0.5, 0.25]
coeff = [1, 0, 0, 0]
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=0_I4B, x=x )
call ok( ALL(SOFTEQ(ans, exact, tol )))
n = 3
x = [-0.5, 0.5, 0.25]
coeff = [0, 0, 0, 1]
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=3_I4B, x=x )
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 = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff )
exact = Chebyshev1GradientEval( n=n, x=x )
call ok( ALL(SOFTEQ(ans, exact, tol )))
end program main
Interface 3
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1GradientEvalSum(n, x, coeff, k) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of 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 Chebyshev1 polynomial of order n at point x
END FUNCTION Chebyshev1GradientEvalSum
END INTERFACE
This example shows the usage of Chebyshev1GradientEvalSum
method.
This routine evaluates the kth gradient of sum of finite series of Chebyshev1 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
n = 3; coeff = [1, 0, 0, 0]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff, k=1 )
exact = Chebyshev1GradientEval( n=0_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 1, 0, 0]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff, k=1 )
exact = Chebyshev1GradientEval( n=1_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 0, 1]
x = 0.25_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff, k=1 )
exact = Chebyshev1GradientEval( n=3_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
n = 3; coeff = [0, 0, 1, 2]
x = -0.5_DFP
ans = Chebyshev1GradientEvalSum( n=n, x=x, coeff=coeff, k = 1 )
exact = Chebyshev1GradientEval( n=2_I4B, x=x ) &
& + 2.0_DFP * Chebyshev1GradientEval( &
& n=3_I4B, x=x )
call ok( SOFTEQ(ans, exact, tol ))
end program main
Interface 4
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1GradientEvalSum(n, x, coeff, k) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of 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 Chebyshev1 polynomial of order n at point x
END FUNCTION Chebyshev1GradientEvalSum
END INTERFACE
See example of Interface 2.