JacobiGaussLobattoQuadrature
This routine returns the Quadrature points and weights.
The Gauss-Lobatto quadrature points consists both as quadrature points.
- The first quadrature point is
- The second quadrature point is
The remaining points are internal to , and they are n-zeros of Jacobi polynomial of order n with respect to the following weight.
Here n is the order of Jacobi polynomial.
INTERFACE
MODULE SUBROUTINE JacobiGaussLobattoQuadrature(n, alpha, beta, pt, wt)
INTEGER(I4B), INTENT(IN) :: n
!! order of Jacobi polynomial
REAL(DFP), INTENT(IN) :: alpha
REAL(DFP), INTENT(IN) :: beta
REAL(DFP), INTENT(OUT) :: pt(:)
!! n+2 quad points indexed from 1 to n+2
REAL(DFP), INTENT(OUT) :: wt(:)
!! n+2 weights, index from 1 to n+2
END SUBROUTINE JacobiGaussLobattoQuadrature
END INTERFACE
Examples
- ️܀ See example
- ↢
- This example shows the usage of
JacobiGaussLobattoQuadrature
method. - This routine returns the quadrature points for Jacobi weights.
In this example (that is, Legendre polynomial) Note that this routine returns n+2 points, the first point is -1 and last point is +1
program main
use easifembase
implicit none
integer( i4b ) :: n
real( dfp ), allocatable :: pt( : ), wt( : )
real( dfp ), parameter :: alpha=0.0, beta=0.0
type(string) :: msg, astr
order=0
n = 0; call callme
See results
Jacobi Gauss Lobatto points, n+2 = 3 alpha=0 beta=0
pt | wt |
---|---|
-1 | 0.33333 |
-4.48639E-17 | 1.3333 |
1 | 0.33333 |
order=1
n = 1; call callme
See results
Jacobi Gauss Lobatto points, n+2 = 3 alpha=0 beta=0
pt | wt |
---|---|
-1 | 0.33333 |
-4.48639E-17 | 1.3333 |
1 | 0.33333 |
order=2.
n = 2; call callme
See results
Jacobi Gauss Lobatto points, n+2 = 4 alpha=0 beta=0
pt | wt |
---|---|
-1 | 0.16667 |
-0.44721 | 0.83333 |
0.44721 | 0.83333 |
1 | 0.16667 |
order=4.
n = 4; call callme
See results
Jacobi Gauss Lobatto points, n+2 = 4 alpha=0 beta=0
pt | wt |
---|---|
-1 | 6.66667E-02 |
-0.76506 | 0.37847 |
-0.28523 | 0.55486 |
0.28523 | 0.55486 |
0.76506 | 0.37847 |
1 | 6.66667E-02 |
contains
subroutine callme
call reallocate( pt, n+2, wt, n+2 )
call JacobiGaussLobattoQuadrature( n=n, alpha=alpha, &
& beta=beta, pt=pt, wt=wt )
msg="Jacobi Gauss Lobatto points, n+2 = " &
& // tostring( n+2 ) // " alpha="//tostring( alpha ) // &
& " beta="//tostring( beta )
call display(msg%chars())
astr = MdEncode( pt .COLCONCAT. wt )
call display( astr%chars(), "" )
end subroutine callme
end program main