JacobiEvalAll
Evaluate Jacobi polynomials from order = 0 to n at single or several points.
Interface 1
INTERFACE
  MODULE PURE FUNCTION JacobiEvalAll(n, alpha, beta, x) RESULT(ans)
    INTEGER(I4B), INTENT(IN) :: n
    REAL(DFP), INTENT(IN) :: alpha
    REAL(DFP), INTENT(IN) :: beta
    REAL(DFP), INTENT(IN) :: x
    REAL(DFP) :: ans(n + 1)
    !! Evaluate Jacobi polynomial of order = 0 to n (total n+1)
    !! at point x
  END FUNCTION JacobiEvalAll
END INTERFACE
Evaluate Jacobi polynomials from order = 0 to n at single points
- N, the highest order polynomial to compute. Note that polynomials 0 through N will be computed.
- alpha, beta are parameters
- x: the point at which the polynomials are to be evaluated.
- ans(1:N+1), the values of the first N+1 Jacobi polynomials at x
Interface 2
INTERFACE
  MODULE PURE FUNCTION JacobiEvalAll(n, alpha, beta, x) RESULT(ans)
    INTEGER(I4B), INTENT(IN) :: n
    REAL(DFP), INTENT(IN) :: alpha
    REAL(DFP), INTENT(IN) :: beta
    REAL(DFP), INTENT(IN) :: x(:)
    REAL(DFP) :: ans(SIZE(x), n + 1)
    !! Evaluate Jacobi polynomial of order = 0 to n (total n+1)
    !! at point x
  END FUNCTION JacobiEvalAll
END INTERFACE
Evaluate Jacobi polynomials from order = 0 to n at several points
- N, the highest order polynomial to compute. Note that polynomials 0 through N will be computed.
- alpha, beta are parameters
- x: the point at which the polynomials are to be evaluated.
- ans(M,1:N+1), the values of the first N+1 Jacobi polynomials at the point
Examples
- ️܀ See example
- ↢
This example shows the usage of JacobiEvalAll method.
This routine evaluates Jacobi polynomial upto order n, for many points.
In this example (that is, Legendre polynomial)
program main
  use easifembase
  implicit none
  integer( i4b ) :: n
  real( dfp ), allocatable :: ans( :, : ), x( : )
  real( dfp ), parameter :: alpha=0.0_DFP, beta=0.0_DFP
  type(string) :: astr
x = [-1.0, 0.0, 1.0]
n = 3; call callme
See results
| P0 | P1 | P2 | P3 | 
|---|---|---|---|
| 1 | -1 | 1 | -1 | 
| 1 | 0 | -0.5 | -0 | 
| 1 | 1 | 1 | 1 | 
contains
subroutine callme
  ans= JacobiEvalAll( n=n, alpha=alpha, beta=beta, x=x )
  astr = MdEncode( ans )
  call display( astr%chars(), "" )
end subroutine callme
end program main
- ️܀ See example
- ↢
This example shows the usage of JacobiEvalAll method.
This routine evaluates Jacobi polynomial upto order n, for many points
In this example (that is, Legendre polynomial)
program main
  use easifembase
  implicit none
  integer( i4b ) :: n
  real( dfp ), allocatable :: ans( : ), x
  real( dfp ), parameter :: alpha=0.0_DFP, beta=0.0_DFP
  type(string) :: astr
x = -1.0_DFP
n = 3; call callme
See results
| P0 | P1 | P2 | P3 | 
|---|---|---|---|
| 1 | -1 | 1 | -1 | 
contains
subroutine callme
  ans= JacobiEvalAll( n=n, alpha=alpha, beta=beta, x=x )
  astr = MdEncode( ans )
  call display( astr%chars(), "" )
end subroutine callme
end program main