UnscaledLobattoEvalAll
Evaluate all UnscaledLobatto polynomials upto order n.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoEvalAll(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans(n + 1)
!! Evaluate UnscaledLobatto polynomial of order = 0 to n (total n+1)
!! at point x
END FUNCTION UnscaledLobattoEvalAll
END INTERFACE
This example shows the usage of UnscaledLobattoEvalAll
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; 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.375 | -0.1875 | -2.34375E-02 | 5.85938E-02 |
contains
subroutine callme
ans= UnscaledLobattoEvalAll( n=n, x=x )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UnscaledLobattoEvalAll(n, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(SIZE(x), n + 1)
!! Evaluate UnscaledLobatto polynomial of order = 0 to n (total n+1)
!! at point x
END FUNCTION UnscaledLobattoEvalAll
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