UltrasphericalQuadrature
This routine returns the Quadrature point of Ultraspherical polynomial.
Here n is the number of quadrature points. Please note it is not the order of Ultraspherical polynomial. The order is decided internally depending upon the quadType
pt
and wt
should be allocated outside, and length should be n.
Interface
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE UltrasphericalQuadrature(n, lambda, 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 Ultraspherical polynomial
!! for quadType = GaussRadauLeft or GaussRadauRight n is order+1
!! for quadType = GaussLobatto, n = order+2
REAL(DFP), INTENT(IN) :: lambda
!! lambda should be greater than -0.5
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 UltrasphericalQuadrature
END INTERFACE
This example shows the usage of UltrasphericalQuadrature
method.
This routine returns the quadrature points for Ultraspherical weights.
By using this subroutine we can get Ultraspherical-Gauss, Ultraspherical-Gauss-Radau, Ultraspherical-Gauss-Lobatto quadrature points
program main
use easifembase
implicit none
integer( i4b ) :: n, quadType
real( dfp ), parameter :: lambda=0.5_DFP
real( dfp ), allocatable :: pt( : ), wt( : )
type(string) :: msg, astr
n = 5; quadType=Gauss; call callme
See results
pt | wt |
---|---|
-0.90618 | 0.23693 |
-0.53847 | 0.47863 |
-1.56541E-16 | 0.56889 |
0.53847 | 0.47863 |
0.90618 | 0.23693 |
n = 5; quadType=GaussRadauLeft; call callme
See results
pt | wt |
---|---|
-1 | 8E-02 |
-0.72048 | 0.44621 |
-0.16718 | 0.62365 |
0.44631 | 0.56271 |
0.88579 | 0.28743 |
n = 5; quadType=GaussRadauRight; call callme
See results
pt | wt |
---|---|
-0.88579 | 0.28743 |
-0.44631 | 0.56271 |
0.16718 | 0.62365 |
0.72048 | 0.44621 |
1 | 8E-02 |
n = 5; quadType=GaussLobatto; call callme
See results
pt | wt |
---|---|
-1 | 0.1 |
-0.65465 | 0.54444 |
-6.41178E-17 | 0.71111 |
0.65465 | 0.54444 |
1 | 0.1 |
contains
subroutine callme
call reallocate( pt, n, wt, n )
call UltrasphericalQuadrature( n=n, pt=pt, wt=wt, &
& quadType=quadType, lambda=lambda )
msg = "| pt | wt |"
call display(msg%chars())
astr = MdEncode( pt .COLCONCAT. wt )
call display( astr%chars(), "" )
end subroutine callme
end program main