LegendreEvalSum
Evaluate the finite sum of Legendre polynomials at point x.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreEvalSum(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 Legendre polynomial of order n at point x
END FUNCTION LegendreEvalSum
END INTERFACE
This example shows the usage of LegendreEvalSum
method.
This routine evaluates the finite sum of the Legendre 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
type(string) :: astr
n = 3
x = [-0.5, 0.5, 0.25]
coeff = [1, 0, 0, 0]
ans = LegendreEvalSum( n=n, x=x, coeff=coeff )
exact = LegendreEval( 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 = LegendreEvalSum( n=n, x=x, coeff=coeff )
exact = LegendreEval( 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 = LegendreEvalSum( n=n, x=x, coeff=coeff )
exact = LegendreEval( n=n, x=x )
call ok( ALL(SOFTEQ(ans, exact, tol )))
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreEvalSum(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 Legendre polynomial of order n at point x
END FUNCTION LegendreEvalSum
END INTERFACE
See above example.