LegendreQuadrature
This routine returns the Quadrature point of Legendre polynomial.
Interface
- ܀ Interface
 - ️܀ See example
 - Example 2
 - ↢
 
INTERFACE
  MODULE SUBROUTINE LegendreQuadrature(n, pt, wt, quadType, onlyInside)
    INTEGER(I4B), INTENT(IN) :: n
    !! number of quadrature points, the order will be computed as follows
    !! for quadType = Gauss, n is same as order of Legendre polynomial
    !! for quadType = GaussRadauLeft or GaussRadauRight n is order+1
    !! for quadType = GaussLobatto, n = order+2
    REAL(DFP), INTENT(OUT) :: pt(n)
    !! n+1 quadrature points from 1 to n+1
    REAL(DFP), OPTIONAL, INTENT(OUT) :: wt(n)
    !! n+1 weights from 1 to n+1
    INTEGER(I4B), INTENT(IN) :: quadType
    !! Gauss
    !! GaussRadauLeft
    !! GaussRadauRight
    !! GaussLobatto
    LOGICAL(LGT), OPTIONAL, INTENT(IN) :: onlyInside
    !! only inside
  END SUBROUTINE LegendreQuadrature
END INTERFACE
nn is the total number of quadrature points.
npoint Gauss quadrature rule has degree of accuracy.npoint Gauss-Radau quadrature rule has degree of accuracy.npoint Gauss-Lobatto quadrature rule has degree of accuracy.
quadTypeType of quadrature, following values are accepted:
- Gauss
 - GaussRadauLeft
 - GaussRadauRight
 - GaussLobatto
 
This example shows the usage of LegendreQuadrature method.
This routine returns the quadrature points for Legendre weights.
By using this subroutine we can get Legendre-Gauss, Legendre-Gauss-Radau, Legendre-Gauss-Lobatto quadrature points
program main
  use easifembase
  implicit none
  integer( i4b ) :: n, quadType
  real( dfp ), allocatable :: pt( : ), wt( : )
  type(string) :: msg, astr
"Legendre-Gauss"
n = 2; quadType=Gauss; call callme
| pt | wt | 
|---|---|
| -0.57735 | 1 | 
| 0.57735 | 1 | 
"Legendre-Radau-Left"
n = 3; quadType=GaussRadauLeft; call callme
| pt | wt | 
|---|---|
| -1 | 0.22222 | 
| -0.2899 | 1.025 | 
| 0.6899 | 0.75281 | 
"Legendre-Radau-Right"
n = 3; quadType=GaussRadauRight; call callme
| pt | wt | 
|---|---|
| -0.6899 | 0.75281 | 
| 0.2899 | 1.025 | 
| 1 | 0.22222 | 
"Legendre-Lobatto"
n = 4; quadType=GaussLobatto; call callme
| pt | wt | 
|---|---|
| -1 | 0.16667 | 
| -0.44721 | 0.83333 | 
| 0.44721 | 0.83333 | 
| 1 | 0.16667 | 
  contains
  subroutine callme
    call reallocate( pt, n, wt, n )
    call LegendreQuadrature( n=n, pt=pt, wt=wt, &
      & quadType=quadType )
    msg = "| pt | wt |"
    call display(msg%chars())
    astr = MdEncode( pt .COLCONCAT. wt )
    call display( astr%chars(), "" )
  end subroutine callme
end program main
This example shows the usage of LegendreQuadrature method.
This routine returns the quadrature points for Legendre polynomials, also we get only inside quadrature points and their weights.
By using this subroutine we can get Legendre-Gauss, Legendre-Gauss-Radau, Legendre-Gauss-Lobatto quadrature points
program main
  use easifembase
  implicit none
  integer( i4b ) :: n, quadType
  real( dfp ), allocatable :: pt( : ), wt( : )
  type(string) :: msg, astr
  logical(lgt), parameter :: onlyInside=.true.
"Legendre-Gauss"
n = 2; quadType=Gauss; call callme
| pt | wt | 
|---|---|
| -0.57735 | 1 | 
| 0.57735 | 1 | 
"Legendre-Radau-Left"
n = 2; quadType=GaussRadauLeft; call callme
| pt | wt | 
|---|---|
| -0.2899 | 1.025 | 
| 0.6899 | 0.75281 | 
"Legendre-Radau-Right"
n = 2; quadType=GaussRadauRight; call callme
| pt | wt | 
|---|---|
| -0.6899 | 0.75281 | 
| 0.2899 | 1.025 | 
"Legendre-Lobatto"
n = 2; quadType=GaussLobatto; call callme
| pt | wt | 
|---|---|
| -0.44721 | 0.83333 | 
| 0.44721 | 0.83333 | 
  contains
  subroutine callme
    call reallocate( pt, n, wt, n )
    call LegendreQuadrature( n=n, pt=pt, wt=wt, &
      & quadType=quadType, onlyInside=onlyInside )
    msg = "| pt | wt |"
    call display(msg%chars())
    astr = MdEncode( pt .COLCONCAT. wt )
    call display( astr%chars(), "" )
  end subroutine callme
end program main