UltrasphericalEvalAll
Evaluate Ultraspherical polynomials from order = 0 to n at single or several points.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEvalAll(n, lambda, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda should be greater than -0.5
REAL(DFP), INTENT(IN) :: x
REAL(DFP) :: ans(n + 1)
!! Evaluate Ultraspherical polynomial of order = 0 to n (total n+1)
!! at point x
END FUNCTION UltrasphericalEvalAll
END INTERFACE
This example shows the usage of UltrasphericalEvalAll
method.
This routine evaluates Ultraspherical polynomial upto order n, for many points
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans( :, : ), x( : )
real( dfp ), parameter :: lambda=0.5_DFP
type(string) :: astr
x = [-1.0, -0.5, 0.0, 0.5, 1.0]
n = 5; call callme
See results
P0 | P1 | P2 | P3 | P4 | P5 |
---|---|---|---|---|---|
1 | -1 | 1 | -1 | 1 | -1 |
1 | -0.5 | -0.125 | 0.4375 | -0.28906 | -8.98438E-02 |
1 | 0 | -0.5 | -0 | 0.375 | 0 |
1 | 0.5 | -0.125 | -0.4375 | -0.28906 | 8.98438E-02 |
1 | 1 | 1 | 1 | 1 | 1 |
contains
subroutine callme
ans= UltrasphericalEvalAll( n=n, x=x, lambda=lambda )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION UltrasphericalEvalAll(n, lambda, x) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: lambda
!! lambda should be greater than -0.5
REAL(DFP), INTENT(IN) :: x(:)
REAL(DFP) :: ans(SIZE(x), n + 1)
!! Evaluate Ultraspherical polynomial of order = 0 to n (total n+1)
!! at point x
END FUNCTION UltrasphericalEvalAll
END INTERFACE
This example shows the usage of UltrasphericalEvalAll
method.
This routine evaluates Ultraspherical polynomial upto order n, at a single point
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: ans( : ), x
real( dfp ), parameter :: lambda=0.5_DFP
type(string) :: astr
x = 0.5_DFP
n = 5; call callme
See results
P0 | P1 | P2 | P3 | P4 | P5 |
---|---|---|---|---|---|
1 | 0.5 | -0.125 | -0.4375 | -0.28906 | 8.98438E-02 |
contains
subroutine callme
ans= UltrasphericalEvalAll( n=n, x=x, lambda=lambda )
astr = MdEncode( ans )
call display( astr%chars(), "" )
end subroutine callme
end program main