Chebyshev1Quadrature
This routine returns the Quadrature point of Chebyshev1 polynomial.
Interface
- ܀ Interface
- ️܀ See example
- Example 2
- ↢
INTERFACE
MODULE SUBROUTINE Chebyshev1Quadrature(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 Chebyshev1 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 Chebyshev1Quadrature
END INTERFACE
By using this subroutine we can get Chebyshev1-Gauss, Chebyshev1-Gauss-Radau, Chebyshev1-Gauss-Lobatto quadrature points
program main
use easifembase
implicit none
integer( i4b ) :: n, quadType
real( dfp ), allocatable :: pt( : ), wt( : )
type(string) :: msg, astr
n = 2; quadType=Gauss; call callme
See results
pt | wt |
---|---|
-0.70711 | 1.5708 |
0.70711 | 1.5708 |
n = 3; quadType=GaussRadauLeft; call callme
See results
pt | wt |
---|---|
-1 | 0.31416 |
-0.30902 | 0.62832 |
0.80902 | 0.62832 |
n = 3; quadType=GaussRadauRight; call callme
See results
pt | wt |
---|---|
-0.80902 | 0.31416 |
0.30902 | 0.62832 |
1 | 0.62832 |
n = 4; quadType=GaussLobatto; call callme
See results
pt | wt |
---|---|
-1 | 0.5236 |
-0.5 | 1.0472 |
0.5 | 1.0472 |
1 | 0.5236 |
contains
subroutine callme
call reallocate( pt, n, wt, n )
call Chebyshev1Quadrature( 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
program main
use easifembase
implicit none
integer( i4b ) :: n, quadType
real( dfp ), allocatable :: pt( : ), wt( : )
type(string) :: msg, astr
logical(lgt), parameter :: onlyInside=.true.
n = 2; quadType=Gauss; call callme
pt | wt |
---|---|
-0.70711 | 1.5708 |
0.70711 | 1.5708 |
n = 2; quadType=GaussRadauLeft; call callme
pt | wt |
---|---|
-0.30902 | 0.62832 |
0.80902 | 0.62832 |
n = 2; quadType=GaussRadauRight; call callme
pt | wt |
---|---|
-0.80902 | 0.31416 |
0.30902 | 0.62832 |
n = 2; quadType=GaussLobatto; call callme
pt | wt |
---|---|
-0.5 | 1.0472 |
0.5 | 1.0472 |
contains
subroutine callme
call reallocate( pt, n, wt, n )
call Chebyshev1Quadrature( 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