LegendreTransform
Discrete Legendre transform.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreTransform(n, coeff, x, w, &
& quadType) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of jacobi polynomial
REAL(DFP), INTENT(IN) :: coeff(0:n)
!! nodal value (at quad points)
REAL(DFP), INTENT(IN) :: x(0:n)
!! quadrature points
REAL(DFP), INTENT(IN) :: w(0:n)
!! weights
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature type, Gauss, GaussLobatto, GaussRadau, GaussRadauLeft
!! GaussRadauRight
REAL(DFP) :: ans(0:n)
!! modal values or coefficients
END FUNCTION LegendreTransform
END INTERFACE
- This example shows the usage of
LegendreTransform
method. - This routine performs the forward Legendre transform.
program main
use easifembase
implicit none
integer( i4b ) :: n
real(dfp), allocatable :: uhat(:), u( : ), pt( : ), &
& wt(:), exact(:)
real( dfp ), parameter :: tol=1.0E-10
type(string) :: astr
integer( i4b ), parameter :: quadType = GaussLobatto
n = 10
call reallocate( pt, n+1, wt, n+1 )
call LegendreQuadrature( n=n+1, pt=pt, wt=wt, &
& quadType=quadType )
u = LegendreEval(n=5_I4B, x=pt)
uhat = LegendreTransform(n=n, coeff=u, x=pt, w=wt, &
& quadType=quadType)
uhat( 6 ) = uhat( 6 ) - 1.0_dFP
CALL ok( SOFTEQ( NORM2(uhat), 0.0_DFP, tol), "test" )
n = 10
call reallocate( pt, n+1, wt, n+1 )
call LegendreQuadrature( n=n+1, &
& pt=pt, wt=wt, quadType=quadType )
u = SIN(4.0_DFP * pi * pt)
uhat = LegendreTransform(n=n, &
& coeff=u, x=pt, w=wt, quadType=quadType)
exact = UltrasphericalTransform(n=n, &
& lambda=0.5_DFP, coeff=u, x=pt, w=wt, &
& quadType=quadType)
call OK(ALL(SOFTEQ(uhat, exact, tol)), "n=10")
n = 20
call reallocate( pt, n+1, wt, n+1 )
call LegendreQuadrature( n=n+1, &
& pt=pt, wt=wt, quadType=quadType )
u = SIN(4.0_DFP * pi * pt)
uhat = LegendreTransform(n=n, &
& coeff=u, x=pt, w=wt, quadType=quadType)
exact = UltrasphericalTransform(n=n, &
& lambda=0.5_DFP, coeff=u, x=pt, w=wt, &
& quadType=quadType)
call OK(ALL(SOFTEQ(uhat, exact, tol)), "n=10")
CALL afile%deallocate()
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION LegendreTransform(n, coeff, x, w, &
& quadType) RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of polynomial
REAL(DFP), INTENT(IN) :: coeff(0:, 1:)
!! nodal value (at quad points)
REAL(DFP), INTENT(IN) :: x(0:n)
!! quadrature points
REAL(DFP), INTENT(IN) :: w(0:n)
!! weights
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature type, Gauss, GaussLobatto, GaussRadau, GaussRadauLeft
!! GaussRadauRight
REAL(DFP) :: ans(0:n, 1:SIZE(coeff, 2))
!! modal values or coefficients for each column of val
END FUNCTION LegendreTransform
END INTERFACE
- This example shows the usage of
LegendreTransform
method. - This routine performs the forward Legendre transform of column data.
program main
use easifembase
implicit none
integer( i4b ) :: n
real(dfp), allocatable :: uhat(:, :), u( :, : ), pt( : ), wt(:), &
& exact( :, : )
real( dfp ), parameter :: tol=1.0E-10
type(string) :: astr
integer( i4b ), parameter :: quadType = GaussLobatto
n = 10
call reallocate( pt, n+1, wt, n+1 )
call LegendreQuadrature( n=n+1, pt=pt, wt=wt, quadType=quadType )
u = SIN(4.0_DFP * pi * pt) .COLCONCAT. SIN(4.0_DFP * pi * pt)
uhat = LegendreTransform(n=n, coeff=u, x=pt, w=wt, quadType=quadType)
exact = UltrasphericalTransform(n=n, lambda=0.5_DFP, coeff=u, x=pt, &
& w=wt, quadType=quadType)
call OK(ALL(SOFTEQ(uhat, exact, tol)), "n=10")
end program main
Interface 3
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE FUNCTION LegendreTransform(n, f, quadType) &
& RESULT(ans)
INTEGER(I4B), INTENT(IN) :: n
!! order of jacobi polynomial
PROCEDURE(iface_1DFunction), POINTER, INTENT(IN) :: f
!! 1D space function
INTEGER(I4B), INTENT(IN) :: quadType
!! Quadrature type, Gauss, GaussLobatto, GaussRadau, GaussRadauLeft
!! GaussRadauRight
REAL(DFP) :: ans(0:n)
!! modal values or coefficients
END FUNCTION LegendreTransform
END INTERFACE
- TODO.