SetAbstractFEParam
Set the parameters necessary for initiating an AbstractFE_
or any of its children.
Interface
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE SetFiniteElementParam
MODULE SUBROUTINE SetAbstractFEParam( &
& param, &
& nsd, &
& elemType, &
& baseContinuity, &
& baseInterpol, &
& ipType, &
& basisType, &
& alpha, &
& beta, &
& lambda, &
& order, &
& anisoOrder, &
& edgeOrder, &
& faceOrder, &
& cellOrder)
TYPE(ParameterList_) :: param
INTEGER(I4B), INTENT(IN) :: nsd
!! Number of spatial dimension
INTEGER(I4B), INTENT(IN) :: elemType
!! Type of finite element
!! Line, Triangle, Quadrangle, Tetrahedron, Prism, Pyramid,
!! Hexahedron
CHARACTER(*), INTENT(IN) :: baseContinuity
!! Continuity or Conformity of basis function.
!! This parameter is used to determine the nodal coordinates of
!! reference element, when xij is not present.
!! If xij is present then this parameter is ignored
!! H1* (default), HDiv, HCurl, DG
CHARACTER(*), INTENT(IN) :: baseInterpol
!! Basis function family used for interpolation.
!! This parameter is used to determine the nodal coordinates of
!! reference element, when xij is not present.
!! If xij is present then this parameter is ignored
!! LagrangeInterpolation, LagrangePolynomial
!! SerendipityInterpolation, SerendipityPolynomial
!! HierarchyInterpolation, HierarchyPolynomial
!! OrthogonalInterpolation, OrthogonalPolynomial
!! HermitInterpolation, HermitPolynomial
INTEGER(I4B), OPTIONAL, INTENT(IN) :: ipType
!! Interpolation point type, It is required when
!! baseInterpol is LagrangePolynomial
INTEGER(I4B), OPTIONAL, INTENT(IN) :: basisType(:)
!! Basis type: Legendre, Lobatto, Ultraspherical,
!! Jacobi, Monomial
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha(:)
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: beta(:)
!! Jacobi parameter
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda(:)
!! Ultraspherical parameters
INTEGER(I4B), OPTIONAL, INTENT(IN) :: order
!! Isotropic Order of finite element
INTEGER(I4B), OPTIONAL, INTENT(IN) :: anisoOrder(:)
!! Anisotropic order, order in x, y, and z directions
INTEGER(I4B), OPTIONAL, INTENT(IN) :: edgeOrder(:)
!! Order of approximation along edges
INTEGER(I4B), OPTIONAL, INTENT(IN) :: faceOrder(:)
!! Order of approximation along face
INTEGER(I4B), OPTIONAL, INTENT(IN) :: cellOrder(:)
!! Order of approximation along cell
END SUBROUTINE SetAbstractFEParam
END INTERFACE SetFiniteElementParam
PROGRAM main
USE easifemBase
USE easifemClasses
IMPLICIT NONE
TYPE(FiniteElement_) :: fe
TYPE(ParameterList_) :: param
INTEGER(I4B) :: nsd, elemType, order, iptype
CALL FPL_Init
CALL param%Initiate()
nsd = 1_I4B
elemType = Line2
order = 1
iptype = Equidistance
CALL SetFiniteElementParam( &
& param=param, &
& nsd=nsd, &
& elemType=elemType, &
& baseContinuity="H1", &
& baseInterpolation="Lagrange", &
& iptype=iptype, &
& basisType=[Monomial], &
& order=order)
CALL fe%Initiate(param)
CALL Display(fe%mdencode(), "")
CALL param%DEALLOCATE()
CALL FPL_Finalize
END PROGRAM main
See results
- Finite Element
- Reference Element
nsd | 1 |
feType | 1 |
elemType | Line2 |
ipType | 1 |
basisType | 0, 0, 0 |
alpha | 0, 0, 0 |
beta | 0, 0, 0 |
lambda | 0, 0, 0 |
dofType | 1, 1, 1, 1 |
transformType | 1 |
baseContinuity | H1 |
baseInterpolion | Lagrange |
refElemDomain | BIUNIT |
isIsotropicOrder | I |
isAnisotropicOrder | I |
isEdgeOrder | I |
isFaceOrder | I |
isCellOrder | I |
edgeOrder | |
faceOrder | |
cellOrder |
Element type | Line2 |
Xidimension | 1 |
NSD | 1 |
tPoints | 2 |
tLines | 1 |
tSurfaces | 0 |
tVolumes | 0 |
BaseContinuity | H1 |
BaseInterpolation | LagrangeInterpolation |
Nodal Coordinates:
x | -1 | 1 |
- PointTopology( 1 ) :
- PointTopology( 2 ) :
Element type | Point1 |
Xidimension | 0 |
Nptrs | 1 |
Element type | Point1 |
Xidimension | 0 |
Nptrs | 2 |
- EdgeTopology( 1 ) :
Element type | Line2 | |
Xidimension | 1 | |
Nptrs | 1 | 2 |
Finite elements on Line
baseContinuity | baseInterpol | baseType | ipType |
---|---|---|---|
H1 | Lagrange | Monomial , Legendre , Lobatto , Jacobi , Ultraspherical | Equidistance , LegendreLobatto , ChebyshevLobatto , UltrasphericalLobatto , JacobiLobatto |
H1 | Orthogonal | Legendre , Lobatto , Jacobi , Ultraspherical | NA |
H1 | Hierarchy | NA | NA |
- 🎉 For
baseType=Jacobi
,alpha
andbeta
should be specified. - 🚀 For
baseType=Ultraspherical
,lambda
should be specified.
Finite elements on Triangle
baseContinuity | baseInterpol | baseType | ipType |
---|---|---|---|
H1 | Lagrange | Monomial , Hierarchy , Orthogonal | Equidistance , LegendreLobatto , ChebyshevLobatto , UltrasphericalLobatto , JacobiLobatto |
H1 | Hierarchy | NA | NA |
H1 | Orthogonal | NA | NA |
Finite elements on Quadrangle
baseContinuity | baseInterpol | baseType | ipType |
---|---|---|---|
H1 | Lagrange | Monomial , Hierarchy , Orthogonal | Equidistance , LegendreLobatto , ChebyshevLobatto , UltrasphericalLobatto , JacobiLobatto |
H1 | Hierarchy | NA | NA |
H1 | Orthogonal | NA | NA |