Chebyshev1Transform
Discrete Chebyshev1 transform.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1Transform(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 Chebyshev1Transform
END INTERFACE
- This example shows the usage of
Chebyshev1Transform
method. - This routine performs the forward Chebyshev1 transform.
program main
use easifembase
implicit none
integer( i4b ) :: n
real(dfp), allocatable :: uhat(:), u( : ), pt( : ), &
& wt(:), exact(:), error(:,:)
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 Chebyshev1Quadrature( n=n+1, pt=pt, wt=wt, &
& quadType=quadType )
u = Chebyshev1Eval(n=5_I4B, x=pt)
uhat = Chebyshev1Transform(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, exact, n+1 )
call Chebyshev1Quadrature( n=n+1, &
& pt=pt, wt=wt, quadType=quadType )
!u = 1.0 / (1.0 + 16.0*pt**2)
u = 1.0*pt + 2.0*pt**2
uhat = Chebyshev1Transform(n=n, &
& coeff=u, x=pt, w=wt, quadType=quadType)
exact(1:3) = 1.0_DFP
call ok( ALL(SOFTEQ(exact, uhat, tol)), "n=10 "//CHAR_LF)
end program main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE PURE FUNCTION Chebyshev1Transform(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 Chebyshev1Transform
END INTERFACE
Interface 3
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE FUNCTION Chebyshev1Transform(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 Chebyshev1Transform
END INTERFACE
- This example shows the usage of
Chebyshev1Transform
method. - This routine performs the forward Chebyshev1 transform.
program main
use easifembase
implicit none
integer( i4b ) :: n
real(dfp), allocatable :: uhat(:), u( : ), pt( : ), &
& wt(:), exact(:), error(:,:)
real( dfp ), parameter :: tol=1.0E-10
type(string) :: astr
integer( i4b ), parameter :: quadType = Gauss
n = 10
call reallocate( pt, n+1, wt, n+1 )
call Chebyshev1Quadrature( n=n+1, pt=pt, wt=wt, &
& quadType=quadType )
u = Chebyshev1Eval(n=5_I4B, x=pt)
uhat = Chebyshev1Transform(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, exact, n+1 )
call Chebyshev1Quadrature( n=n+1, &
& pt=pt, wt=wt, quadType=quadType )
!u = 1.0 / (1.0 + 16.0*pt**2)
u = 1.0*pt + 2.0*pt**2
uhat = Chebyshev1Transform(n=n, &
& coeff=u, x=pt, w=wt, quadType=quadType)
exact(1:3) = 1.0_DFP
call ok( ALL(SOFTEQ(exact, uhat, tol)), "n=10 "//CHAR_LF)
end program main