Chebyshev1InvTransform
Discrete Inverse Chebyshev1 transform.
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
  MODULE PURE FUNCTION Chebyshev1InvTransform(n, coeff, x) &
        & RESULT(ans)
    INTEGER(I4B), INTENT(IN) :: n
    !! order of Jacobi polynomial
    REAL(DFP), INTENT(IN) :: coeff(0:n)
    !! n+1  coefficient (modal values)
    REAL(DFP), INTENT(IN) :: x
    !! x point in physical space
    REAL(DFP) :: ans
    !! value in physical space
  END FUNCTION Chebyshev1InvTransform
END INTERFACE
- This example shows the usage of Chebyshev1InvTransformmethod.
- This routine performs the inverse Chebyshev1 transform at a single point.
program main
  use easifembase
  implicit none
  integer( i4b ) :: n
  real(dfp), allocatable :: uhat(:), u( : ), pt( : ), wt(:), error(:,:)
  real( dfp ), parameter :: tol=1.0E-10
  type(string) :: astr
  integer( i4b ), parameter :: quadType = GaussLobatto
  real( dfp ) :: x, ans, exact
  call reallocate(error, 30, 2)
  do n = 1, size(error,1)
    call reallocate( pt, n+1, wt, n+1 )
    call Chebyshev1Quadrature( n=n+1, pt=pt, wt=wt, quadType=quadType )
    u = 1.0 / (1.0 + 16.0 * pt**2)
    uhat = Chebyshev1Transform(n=n, coeff=u, x=pt, w=wt, quadType=quadType)
    x = 0.1
    exact = 1.0 / (1.0 + 16.0 * x**2)
    ans = Chebyshev1InvTransform(n=n, x=x, coeff=uhat)
    error(n,1) = n
    error(n,2) = ABS(exact-ans)
  end do
  call display(MdEncode(error), "error"//CHAR_LF )
end program main
| n | error | 
|---|---|
| 1 | 0.80325 | 
| 2 | 0.12852 | 
| 3 | 0.61689 | 
| 4 | 0.11195 | 
| 5 | 0.3907 | 
| 10 | 4.75496E-02 | 
| 15 | 2.77239E-03 | 
| 20 | 4.28428E-03 | 
| 25 | 2.75076E-03 | 
| 30 | 5.40866E-05 | 
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
  MODULE PURE FUNCTION Chebyshev1InvTransform(n, coeff, x) &
        & RESULT(ans)
    INTEGER(I4B), INTENT(IN) :: n
    !! order of Jacobi polynomial
    REAL(DFP), INTENT(IN) :: coeff(0:n)
    !! n+1  coefficient (modal values)
    REAL(DFP), INTENT(IN) :: x(:)
    !! x point in physical space
    REAL(DFP) :: ans(SIZE(x))
    !! value in physical space
  END FUNCTION Chebyshev1InvTransform
END INTERFACE
- This example shows the usage of Chebyshev1InvTransformmethod.
- This routine performs the inverse Chebyshev1 transform at several points.
program main
  use easifembase
  use easifemclasses
  implicit none
  integer( i4b ) :: ii, n
  real(dfp), allocatable :: uhat(:), u( : ), pt( : ), wt(:), &
    & x(:), ans(:), exact(:)
  real( dfp ), parameter :: tol=1.0E-10
  type(string) :: astr
  integer( i4b ), parameter :: quadType = GaussLobatto
  type(PlPlot_) :: aplot
  character( len=* ), parameter :: fname = "./results/test24"
  CALL aplot%Initiate()
  CALL aplot%Set( &
    & device="svg", &
    & filename=fname//"-%n.svg")
  CALL aplot%figure()
  x = linspace(-1.0_DFP, 1.0_DFP, 101_I4B)
  ii = 0
  do n = 5,25,5
    ii = ii + 1
    CALL aplot%subplot(5,1,ii)
    CALL aplot%setXYLim([-1.0_DFP, 1.0_DFP], [ -2.0_DFP, 2.0_DFP])
    CALL aplot%setTicks()
    call reallocate( pt, n+1, wt, n+1 )
    call Chebyshev1Quadrature( n=n+1, pt=pt, wt=wt, quadType=quadType )
    u = 1.0/(1.0 + 16.0 * pt**2)
    uhat = Chebyshev1Transform(n=n, coeff=u, x=pt, w=wt, quadType=quadType)
    exact = 1.0/(1.0 + 16.0 * x**2)
    ans = Chebyshev1InvTransform(n=n, x=x, coeff=uhat)
    CALL aplot%plot2D( x=x,y=ans)
    CALL aplot%plot2D( x=x,y=exact, pointType=PS_DOT, lineWidth=0.0_DFP )
    CALL aplot%setLabels("x","u(x)","n="//tostring(n))
  end do
  CALL aplot%Show()
  CALL aplot%deallocate()
end program main