JacobiEvalSum
Evaluate the finite sum of Jacobi polynomials at point x.
Interface 1
INTERFACE
MODULE PURE FUNCTION JacobiEvalSum(n, alpha, beta, x, coeff) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: alpha
!! alpha of Jacobi polynomial
REAL(DFP), INTENT(IN) :: beta
!! beta of Jacobi 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 Jacobi polynomial of order n at point x
END FUNCTION JacobiEvalSum
END INTERFACE
Interface 2
INTERFACE
MODULE PURE FUNCTION JacobiEvalSum(n, alpha, beta, x, coeff) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: alpha
!! alpha of Jacobi polynomial
REAL(DFP), INTENT(IN) :: beta
!! beta of Jacobi 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 Jacobi polynomial of order n at point x
END FUNCTION JacobiEvalSum
END INTERFACE
Example
- ️܀ See example
- ↢
This example shows the usage of JacobiEvalSum
method.
This routine evaluates the finite sum of Jacobi polynomial of order upto n at a given point.
In this example (that is, proportional to Legendre polynomial)
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ) :: ans, x, exact
real(dfp), allocatable :: coeff(:)
real( dfp ), parameter :: alpha=0.0_DFP, beta=0.0_DFP, tol=1.0E-10
type(string) :: astr
n = 3
coeff = [1.0, 0.0, 0.0, 0.0]
x = 0.5_DFP;
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=0_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( SOFTEQ(ans, exact, tol ))
n = 3
coeff = [0.0, 1.0, 0.0, 0.0]
x = 0.5_DFP;
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=1_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( SOFTEQ(ans, exact, tol ))
n = 3
coeff = [0.0, 0.0, 0.0, 1.0]
x = 0.5_DFP;
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=3_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( SOFTEQ(ans, exact, tol ))
n = 3
coeff = [1.0, 2.0, 3.0, 4.0]
x = 0.5_DFP;
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = &
& 1.0*JacobiEval(n=0_I4B, x=x, alpha=alpha,beta=beta) &
&+ 2.0*JacobiEval(n=1_I4B, x=x, alpha=alpha,beta=beta) &
&+ 3.0*JacobiEval(n=2_I4B, x=x, alpha=alpha,beta=beta) &
&+ 4.0*JacobiEval(n=3_I4B, x=x, alpha=alpha,beta=beta)
call ok( SOFTEQ(ans, exact, tol ))
end program main
- ️܀ See example
- ↢
This example shows the usage of JacobiEvalSum
method.
This routine evaluates the finite sum of Jacobi polynomial of order upto n at several points.
In this example (that is, proportional to Legendre polynomial)
program main
use easifembase
implicit none
integer( i4b ) :: n
real(dfp), allocatable :: coeff(:), ans(:), x(:), exact(:)
real( dfp ), parameter :: alpha=0.0_DFP, beta=0.0_DFP, tol=1.0E-10
type(string) :: astr
n = 3
coeff = [1.0, 0.0, 0.0, 0.0]
x = [0.5_DFP, -0.5_DFP];
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=0_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( ALL(SOFTEQ(ans, exact, tol )))
n = 3
coeff = [0.0, 1.0, 0.0, 0.0]
x = [0.5_DFP, -0.5_DFP];
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=1_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( ALL(SOFTEQ(ans, exact, tol )))
n = 3
coeff = [0.0, 0.0, 0.0, 1.0]
x = [0.5_DFP, -0.5_DFP];
ans = JacobiEvalSum(n=n, x=x, alpha=alpha, &
& beta=beta, coeff=coeff)
exact = JacobiEval(n=3_I4B, x=x, alpha=alpha,&
& beta=beta)
call ok( ALL(SOFTEQ(ans, exact, tol )))
end program main