LegendreGradientEval
Evaluate gradient of Legendre polynomials.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreGradientEval(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans
END FUNCTION LegendreGradientEval
END INTERFACE
This example shows the usage of LegendreGradientEval
method.
This routine evaluates gradient of Legendre polynomial of order n, at single point.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ) :: ans, x, exact
real( dfp ), parameter :: tol=1.0E-10
n = 3
x = -1.0_DFP; call callme
exact = 0.5_DFP*(3*5.0 * x**2 - 3.0)
call ok( SOFTEQ(ans, exact, tol ))
x = 0.0_DFP; call callme
exact = 0.5_DFP*(3*5.0 * x**2 - 3.0)
call ok( SOFTEQ(ans, exact, tol ))
x = +1.0_DFP; call callme
exact = 0.5_DFP*(3*5.0 * x**2 - 3.0)
call ok( SOFTEQ(ans, exact, tol ))
contains
subroutine callme
ans= LegendreGradientEval( n=n, x=x )
end subroutine callme
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreGradientEval(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(1:SIZE(x))
END FUNCTION LegendreGradientEval
END INTERFACE
See above example.