JacobiGradientEval
Evaluate gradient of Jacobi polynomials.
Interface 1
INTERFACE
MODULE PURE FUNCTION JacobiGradientEval(n, alpha, beta, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of Jacobi polynomial
REAL(DFP), INTENT(IN) :: alpha
!! alpha > -1.0
REAL(DFP), INTENT(IN) :: beta
!! beta > -1.0
REAL(DFP), INTENT(IN) :: x
!! point
REAL(DFP) :: ans
!! Derivative of Jacobi polynomial of order n at point x
END FUNCTION JacobiGradientEval
END INTERFACE
Interface 2
INTERFACE
MODULE PURE FUNCTION JacobiGradientEval(n, alpha, beta, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: alpha
REAL(DFP), INTENT(IN) :: beta
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(SIZE(x))
!! Derivative of Jacobi polynomial of order n at x
END FUNCTION JacobiGradientEval
END INTERFACE
Examples
- ️܀ See example
- ↢
This example shows the usage of JacobiGradientEval
method.
In this example (that is, proportional to Legendre polynomial)
program main
use easifembase
use easifemclasses
implicit none
integer( i4b ) :: n, ii
real( dfp ) :: exact, ans, x
real( dfp ), parameter :: alpha=0.0_DFP, beta=0.0_DFP, tol=1.0E-10
type(string) :: astr
!!
n=5
x = 0.25
exact = LegendreGradientEval(n=n, x=x)
ans = JacobiGradientEval(n=n, alpha=alpha, beta=beta, x=x)
call OK( SOFTEQ( exact, ans, tol), "")
!!
n=10
x = 0.25
exact = LegendreGradientEval(n=n, x=x)
ans = JacobiGradientEval(n=n, alpha=alpha, beta=beta, x=x)
call OK( SOFTEQ( exact, ans, tol), "")
!!
end program main