LobattoGradientEvalAll
Evaluate gradient of Lobatto polynomials upto order n.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LobattoGradientEvalAll(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans(1:n + 1)
END FUNCTION LobattoGradientEvalAll
END INTERFACE
This example shows the usage of LobattoGradientEvalAll
method.
This routine evaluates the first derivative of 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, 0.0, 1.0]
n = 5; call callme
contains
subroutine callme
ans= LobattoGradientEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main
-0.5 | 0.5 | -1.2247 | 1.5811 | -1.8708 | 2.1213 |
-0.5 | 0.5 | 0 | -0.79057 | -0 | 0.7955 |
-0.5 | 0.5 | 1.2247 | 1.5811 | 1.8708 | 2.1213 |
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LobattoGradientEvalAll(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(1:SIZE(x), 1:n + 1)
END FUNCTION LobattoGradientEvalAll
END INTERFACE
This example shows the usage of LobattoGradientEvalAll
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
contains
subroutine callme
ans= LobattoGradientEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main
-0.5 | 0.5 | 1.2247 | 1.5811 | 1.8708 | 2.1213 |