UnscaledLobattoGradientEvalAll
Evaluate gradient of UnscaledLobatto polynomials upto order n.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoGradientEvalAll(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans(1:n + 1)
END FUNCTION UnscaledLobattoGradientEvalAll
END INTERFACE
This example shows the usage of UnscaledLobattoGradientEvalAll
method.
This routine evaluates the first derivative of UnscaledLobatto 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= UnscaledLobattoGradientEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main
-0.5 | 0.5 | -1 | 1 | -1 | 1 |
-0.5 | 0.5 | 0 | -0.5 | -0 | 0.375 |
-0.5 | 0.5 | 1 | 1 | 1 | 1 |
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoGradientEvalAll(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 UnscaledLobattoGradientEvalAll
END INTERFACE
This example shows the usage of UnscaledLobattoGradientEvalAll
method.
This routine evaluates UnscaledLobatto polynomial upto order n, at a given point.
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans( : ), x
type(string) :: astr
x = 1.0_DFP
n = 5
ans= UnscaledLobattoGradientEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end program main
-0.5 | 0.5 | 1 | 1 | 1 | 1 |