UnscaledLobattoGradientEval
Evaluate gradient of UnscaledLobatto polynomials.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoGradientEval(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans
END FUNCTION UnscaledLobattoGradientEval
END INTERFACE
This example shows the usage of UnscaledLobattoGradientEval
method.
This routine evaluates gradient of UnscaledLobatto 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-8
n = 5
x = -1.0_DFP; call callme
exact =(1.0/9.0/8.0)*(315.0*x**4 - 270.0*x**2 + 27.0)
call ok( SOFTEQ(ans, exact, tol ))
x = 0.0_DFP; call callme
exact =(1.0/9.0/8.0)*(315.0*x**4 - 270.0*x**2 + 27.0)
call ok( SOFTEQ(ans, exact, tol ))
x = +1.0_DFP; call callme
exact =(1.0/9.0/8.0)*(315.0*x**4 - 270.0*x**2 + 27.0)
call ok( SOFTEQ(ans, exact, tol ))
contains
subroutine callme
ans= UnscaledLobattoGradientEval( n=n, x=x )
end subroutine callme
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoGradientEval(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(1:SIZE(x))
END FUNCTION UnscaledLobattoGradientEval
END INTERFACE
This example shows the usage of LobattoEvalAll
method.
This routine evaluates Lobatto polynomial upto order n, for many points.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans( : ), x
type(string) :: astr
x = -1.0_DFP
n = 5; call callme
l0 | l1 | l2 | l3 | l4 | l5 |
---|---|---|---|---|---|
1 | 0 | 0 | 0 | 0 | 0 |
x = 1.0_DFP
n = 5; call callme
l0 | l1 | l2 | l3 | l4 | l5 |
---|---|---|---|---|---|
0 | 1 | 0 | 0 | 0 | 0 |
x = 0.5_DFP
n = 5; call callme
l0 | l1 | l2 | l3 | l4 | l5 |
---|---|---|---|---|---|
0.25 | 0.75 | -0.45928 | -0.29646 | -4.38475E-02 | 0.1243 |
contains
subroutine callme
ans= LobattoEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main