Methods
Chebyshev polynomial of first kind is defined. The nth order Chebyshev polynomial is denoted by . Chebyshev polynomials are orthogonal polynomials.
It is a child of [[AbstractOrthopol1D_]].
Structure
TYPE, EXTENDS( AbstractOrthoPol1D_ ) :: ChebyshevFirst1D_
CONTAINS
ConstructorMethods
Constructor function
We can create an instance of by calling ChebyshevFirst1D()
function. Examples are included in
- [[ChebyshevFirst1D_test_1]]
- [[ChebyshevFirst1D_test_2]]
Interface is given below:
MODULE FUNCTION ChebyshevFirst1D( varname, n, isMonic, &
& isOrthonormal ) RESULT( ans )
CHARACTER( LEN = * ), INTENT( IN ) :: varname
!! variable name
INTEGER( I4B ), INTENT( IN ) :: n
!! order of chebyshev polynomial
LOGICAL( LGT ), OPTIONAL, INTENT( IN ) :: isMonic
!! Default is .FALSE., if true then leading coeff of poly is 1
LOGICAL( LGT ), OPTIONAL, INTENT( IN ) :: isOrthonormal
!! Default is .FALSE., if true then the polynomials are orthonormal
TYPE( ChebyshevFirst1D_ ) :: ans
!! Chebyshev polynomial of first kind
END FUNCTION ChebyshevFirst1D
ChebyshevFirst1D_Pointer
We can create an instance of new pointer by ChebyshevFirst1D_Pointer()
function.
The interface is given below:
FUNCTION ChebyshevFirst1D_Pointer( varname, n, &
& isMonic, isOrthonormal ) RESULT( ans )
CHARACTER( LEN = * ), INTENT( IN ) :: varname
!! variable name
INTEGER( I4B ), INTENT( IN ) :: n
!! order of chebyshev polynomial
LOGICAL( LGT ), OPTIONAL, INTENT( IN ) :: isMonic
!! Default is .FALSE., if true then leading coeff of poly is 1
LOGICAL( LGT ), OPTIONAL, INTENT( IN ) :: isOrthonormal
!! Default is .FALSE., if true then the polynomials are orthonormal
CLASS( ChebyshevFirst1D_ ), POINTER :: ans
!! Chebyshev polynomial of first kind
END FUNCTION ChebyshevFirst1D_Pointer
Deallocate
You can deallocate (or destroy) the object by calling Deallocate()
subroutine.
CALL obj%Deallocate()
GetMethods
GetStringForUID
Weight
This function evaluate the weight of the Chebyshev polynomial at a given point x.
GetRecurrenceCoeff
GetCoeffScale
Zeros
GaussQuadrature
GaussRadauQuadratur
GaussLobattoQuadrature
Examples
- 🌻 Example 1
- 🌻 Example 2
- 🌻 Example 3
- 🌻 Example 4
- 🌻 Example 5
- 🔥 Clear
ChebyshevFirst1D example 1
This example shows the usage of [[ChebyshevFirst1D_]] class.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
type(String) :: astr
integer(i4b) :: ii
n=1
obj = ChebyshevFirst1D(varname="x", n=1)
call Display( "T(n=1) := " )
call obj%Display( "=>" )
>result
J(n=1, alpha=0.0, beta=0.0) :=
=>J_1( x ) = ( x-.000 )* 1.000* J_0
=>J_0( x ) = 1
n=2
call blanklines( nol=5 )
obj = ChebyshevFirst1D(varname="x", n=2)
call Display( "J(n=2, alpha=0.0, beta=0.0) := " )
call obj%Display( "=>" )
>result
J(n=2, alpha=0.0, beta=0.0) :=
=>J_2( x ) = ( x-.000 )* 1.500* J_1( x ) - .333* 1.500* J_0( x )
=>J_1( x ) = ( x-.000 )* 1.000* J_0
=>J_0( x ) = 1
n=3
call blanklines( nol=5 )
obj = ChebyshevFirst1D(varname="x", n=3)
call Display( "J(n=3, alpha=0.0, beta=0.0) := " )
call obj%Display( "=>" )
>result
J(n=3, alpha=0.0, beta=0.0) :=
=>J_3( x ) = ( x-.000 )* 1.667* J_2( x ) - .267* 2.500* J_1( x )
=>J_2( x ) = ( x-.000 )* 1.500* J_1( x ) - .333* 1.500* J_0( x )
=>J_1( x ) = ( x-.000 )* 1.000* J_0
=>J_0( x ) = 1
n=4
call blanklines( nol=5 )
obj = Jacobi1D(varname="x", n=4, alpha=0.0_DFP, beta=0.0_DFP)
call Display( "J(n=4, alpha=0.0, beta=0.0) := " )
call obj%Display( "=>" )
>result
J(n=4, alpha=0.0, beta=0.0) :=
=>J_4( x ) = ( x-.000 )* 1.750* J_3( x ) - .257* 2.917* J_2( x )
=>J_3( x ) = ( x-.000 )* 1.667* J_2( x ) - .267* 2.500* J_1( x )
=>J_2( x ) = ( x-.000 )* 1.500* J_1( x ) - .333* 1.500* J_0( x )
=>J_1( x ) = ( x-.000 )* 1.000* J_0
=>J_0( x ) = 1
END PROGRAM main
ChebyshevFirst1D example 2
This example shows the usage of [[ChebyshevFirst1D_]] class. It checks the memory leak.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
type(String) :: astr
integer(i4b) :: ii
Let us generate 100 Chebyshev1 polynomials
DO ii = 1, 1000000
obj=ChebyshevFirst1D(varname="x", n=100)
call obj%Deallocate()
END DO
Read( *, * ) ii
END PROGRAM main
ChebyshevFirst1D example 3
This example shows the usage of [[ChebyshevFirst1D_]] class. We test Eval
in this routine
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
type(String) :: astr
integer(i4b) :: ii
real( dfp ), allocatable :: x(:), yexact( : ), y( : )
real( dfp ), parameter :: tol=1.0E-10
n=1
x = linspace(-1.0_DFP, 1.0_DFP, 11_I4B)
yexact = x
obj = ChebyshevFirst1D(varname="x", n=1)
y = obj%Eval( x )
call OK( ALL(SOFTEQ(y, yexact, tol)), "test-1:" )
n=2
yexact = 2.0*x**2 - 1.0
obj = ChebyshevFirst1D(varname="x", n=2)
y = obj%Eval( x )
call OK( ALL(SOFTEQ(y, yexact, tol)), "test-2:" )
n=3
yexact = 4.0*x**3 - 3.0*x
obj = ChebyshevFirst1D(varname="x", n=3)
y = obj%Eval( x )
call OK( ALL(SOFTEQ(y, yexact, tol)), "test-3:" )
n=4
yexact = 8.0*x**4 - 8.0*x**2 + 1.0
obj = ChebyshevFirst1D(varname="x", n=4)
y = obj%Eval( x )
call OK( ALL(SOFTEQ(y, yexact, tol)), "test-4:" )
END PROGRAM main
ChebyshevFirst1D example 4
This example shows the usage of [[ChebyshevFirst1D_]] class. We test EvalGradient
in this routine. the argument of EvalGradient
is scalar in this routine.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
type(String) :: astr
integer(i4b) :: ii, n
real( dfp ) :: x, yexact, y
real( dfp ), parameter :: tol=1.0E-10
n=1
x = RAND()
yexact = 1
n=1
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK( SOFTEQ( y, yexact, tol ), "test-1:" )
n=2
x = RAND()
yexact = 4.0*x
n=2
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK( SOFTEQ( y, yexact, tol ), "test-2:" )
n=3
x = RAND()
yexact = 12.0*x**2 - 3.0
n=3
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK( SOFTEQ( y, yexact, tol ), "test-3:" )
n=4
x = RAND()
yexact = 32.0*x**3 - 16.0*x
n=4
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK( SOFTEQ( y, yexact, tol ), "test-4:" )
END PROGRAM main
ChebyshevFirst1D example 5
This example shows the usage of [[ChebyshevFirst1D_]] class. We test EvalGradient
in this routine. the argument of EvalGradient
is vector in this routine.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
type(String) :: astr
integer(i4b) :: ii, n
real( dfp ), allocatable :: x(:), yexact(:), y(:)
real( dfp ) :: tol=1.0E-10
n=1
x = linspace(0.0_DFP, 1.0_DFP, 11_I4B)
yexact = ones(11, 1)
n=1
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK( ALL( y .APPROXEQ. yexact), "test-1:" )
n=2
yexact = 4.0*x
n=2
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK(ALL(y .APPROXEQ. yexact), "test-2:" )
n=3
yexact = 12.0*x**2 - 3.0
n=3
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK(ALL( SOFTEQ( y, yexact, tol=tol )), "test-3:" )
n=4
yexact = 32.0*x**3 - 16.0*x
n=4
obj = ChebyshevFirst1D(varname="x", n=n)
y = obj%EvalGradient( x )
call OK(ALL( SOFTEQ( y, yexact, tol=tol )), "test-4:" )
END PROGRAM main
- 🌻 Example 6
- 🌻 Example 7
- 🌻 Example 8
- 🌻 Example 9
- 🔥 Clear
ChebyshevFirst1D example 6
This example shows the usage of [[ChebyshevFirst1D_]] class. We test Zeros
function in this routine, which returns the zeros of Chebyshev1 polynomial.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
real( dfp ), allocatable :: x( : )
integer( i4b ) :: n
n=1
n = 1
obj = ChebyshevFirst1D(varname="x", n=n)
x = obj%zeros()
call display( x, "zeros for n="//tostring(n), orient="ROW" )
n=2
n = 2
obj = ChebyshevFirst1D(varname="x", n=n)
x = obj%zeros()
call display( x, "zeros for n="//tostring(n), orient="ROW" )
n=3
n = 3
obj = ChebyshevFirst1D(varname="x", n=n)
x = obj%zeros()
call display( x, "zeros for n="//tostring(n), orient="ROW" )
n=4
n = 4
obj = ChebyshevFirst1D(varname="x", n=n)
x = obj%zeros()
call display( x, "zeros for n="//tostring(n), orient="ROW" )
END PROGRAM main
ChebyshevFirst1D example 7
This example shows the usage of [[ChebyshevFirst1D_]] class. We test GaussQuadrature
function in this routine, which returns the GaussQuadrature points for Chebyshev1 polynomial.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
real( dfp ), allocatable :: x( :, : )
integer( i4b ) :: n
n=1
n = 1
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=2
n = 2
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=3
n = 3
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=4
n = 4
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
END PROGRAM main
ChebyshevFirst1D example 8
This example shows the usage of [[ChebyshevFirst1D_]] class. We test GaussRadauQuadrature
function in this routine, which returns the GaussRadauQuadrature points for Chebyshev1 polynomial.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
real( dfp ), allocatable :: x( :, : )
integer( i4b ) :: n
real( dfp ), parameter :: a = -1.0_DFP
n=1
n = 1
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussRadauQuadrature(a)
call display( x, "pt | wt for n="//tostring(n) )
n=2
n = 2
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussRadauQuadrature(a)
call display( x, "pt | wt for n="//tostring(n) )
n=3
n = 3
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussRadauQuadrature(a)
call display( x, "pt | wt for n="//tostring(n) )
n=4
n = 4
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussRadauQuadrature(a)
call display( x, "pt | wt for n="//tostring(n) )
END PROGRAM main
ChebyshevFirst1D example 9
This example shows the usage of [[ChebyshevFirst1D_]] class. We test GaussLobattoQuadrature
function in this routine, which returns the Gauss-Lobatto Quadrature points for Chebyshev1 polynomial.
Modules and classes
- [[ChebyshevFirst1D_]]
Usage
PROGRAM main
use easifemBase
use easifemClasses
implicit none
type(ChebyshevFirst1D_) :: obj
real( dfp ), allocatable :: x( :, : )
integer( i4b ) :: n
n=1
n = 1
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussLobattoQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=2
n = 2
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussLobattoQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=3
n = 3
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussLobattoQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
n=4
n = 4
obj=ChebyshevFirst1D(varname="x", n=n)
x = obj%GaussLobattoQuadrature()
call display( x, "pt | wt for n="//tostring(n) )
END PROGRAM main