JacobiGradientEvalAll
Evaluate gradient of Jacobi polynomials.
Interface 1
INTERFACE
MODULE PURE FUNCTION JacobiGradientEvalAll(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(n + 1)
!! Derivative of Jacobi polynomial of order n at point x
END FUNCTION JacobiGradientEvalAll
END INTERFACE
Interface 2
INTERFACE
MODULE PURE FUNCTION JacobiGradientEvalAll(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), n + 1)
!! Derivative of Jacobi polynomial of order n at x
END FUNCTION JacobiGradientEvalAll
END INTERFACE
Examples
- ️܀ See example
- ↢
This example shows the usage of JacobiGradientEvalAll
method.
In this example (that is, proportional to Legendre polynomial)
Module and classes
program main
use easifembase
use easifemclasses
implicit none
integer( i4b ) :: n, ii
real( dfp ), allocatable :: exact(:), ans(:)
real( dfp ) :: 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 = LegendreGradientEvalAll(n=n, x=x)
ans = JacobiGradientEvalAll(n=n, alpha=alpha, beta=beta, x=x)
call OK( ALL(SOFTEQ( exact, ans, tol)), "")
!!
n=10
x = 0.25
exact = LegendreGradientEvalAll(n=n, x=x)
ans = JacobiGradientEvalAll(n=n, alpha=alpha, beta=beta, x=x)
call OK( ALL(SOFTEQ( exact, ans, tol)), "")
!!
end program main
- ️܀ See example
- ↢
This example shows the usage of JacobiGradientEvalAll
method.
In this example (that is, proportional to Legendre polynomial)
program main
use easifembase
use easifemclasses
implicit none
integer( i4b ) :: n, ii
real( dfp ), allocatable :: 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.1, 0.0, 0.1]
exact = LegendreGradientEvalAll(n=n, x=x)
ans = JacobiGradientEvalAll(n=n, alpha=alpha, beta=beta, x=x)
call OK( ALL(SOFTEQ( exact, ans, tol)), "")
!!
n=10
x = [-0.1, 0.0, 0.1]
exact = LegendreGradientEvalAll(n=n, x=x)
ans = JacobiGradientEvalAll(n=n, alpha=alpha, beta=beta, x=x)
call OK( ALL(SOFTEQ( exact, ans, tol)), "")
!!
end program main