Skip to main content

Initiate

This method initiates an instance of FEDOF. There are several ways to initiate an instance of FEDOF.

The Initiate method has four different implementations:

  1. Initiate1 - Initializes with homogeneous order for all elements
  2. Initiate2 - Initializes with inhomogeneous orders specified per element
  3. Initiate3 - Initializes from a parameter list
  4. Initiate4 - Initializes from an order vector defined for global elements

Interface 1

This method is for Homogeneous order, that is, all elements in the mesh have the same order.

INTERFACE
MODULE SUBROUTINE Initiate(obj, order, mesh, baseContinuity, &
baseInterpolation, ipType, basisType, alpha, beta, lambda)
CLASS(FEDOF_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: order
!! homogeneous value of order
CLASS(AbstractMesh_), TARGET, INTENT(IN) :: mesh
!! cell mesh
CHARACTER(*), INTENT(IN) :: baseContinuity
!! continuity of basis (regularity or conformity)
CHARACTER(*), INTENT(IN) :: baseInterpolation
!! basis function used for interpolation
INTEGER(I4B), OPTIONAL, INTENT(IN) :: ipType
!! interpolation point type
!! used when baseInterpolation is Lagrange
INTEGER(I4B), OPTIONAL, INTENT(IN) :: basisType(:)
!! type of basis function used for
!! constructing the Lagrange polynomial
!! Used when baseInterpolation is Lagrange
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha(:)
!! alpha parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: beta(:)
!! beta parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda(:)
!! lambda parameter for Ultraspherical parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Ultraspherical
END SUBROUTINE Initiate
END INTERFACE

Description

This subroutine configures a FEDOF object by setting up the basis functions, continuity, and interpolation properties based on a given mesh and order.

Parameters

ParameterTypeIntentDescription
objCLASS(FEDOF_)INOUTThe FEDOF object to be initialized.
orderINTEGER(I4B)INHomogeneous polynomial order for the basis functions.
meshCLASS(AbstractMesh_), TARGETINThe computational mesh defining the spatial discretization.
baseContinuityCHARACTER(*)INSpecifies the continuity of basis functions (regularity or conformity).
baseInterpolationCHARACTER(*)INType of basis function used for interpolation.

Optional Parameters

ParameterTypeIntentOptionalDescription
ipTypeINTEGER(I4B)INYesInterpolation point type. Used when baseInterpolation is set to Lagrange.
basisTypeINTEGER(I4B)INYesArray specifying types of basis functions used for constructing Lagrange polynomials. Used when baseInterpolation is Lagrange.
alphaREAL(DFP)INYesAlpha parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
betaREAL(DFP)INYesBeta parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
lambdaREAL(DFP)INYesLambda parameters for Ultraspherical polynomials. Required when baseInterpolation is Lagrange and basisType is Ultraspherical.

Usage Example

CALL myFEDOF%Initiate(order=2, mesh=myMesh, baseContinuity="H1", &
baseInterpolation="Lagrange", ipType=GAUSS_LOBATTO)

Interface 2

  • Here order represents the order of each cell element.
  • order is a vector of integers, the length of order must be equal to the number of elements in the mesh.
  • order(i) is the order of local element i.
INTERFACE
MODULE SUBROUTINE Initiate(obj, order, mesh, baseContinuity, &
baseInterpolation, ipType, basisType, alpha, lambda, beta, islocal)
CLASS(FEDOF_), INTENT(INOUT) :: obj
!! Finite degree of freedom object
INTEGER(I4B), INTENT(IN) :: order(:)
!! Inhomogeneous value of order
!! This is order of each cell element
!! see the note on islocal
CLASS(AbstractMesh_), TARGET, INTENT(IN) :: mesh
!! cell mesh
CHARACTER(*), INTENT(IN) :: baseContinuity
!! continuity of basis (regularity)
CHARACTER(*), INTENT(IN) :: baseInterpolation
!! basis function used for interpolation
INTEGER(I4B), OPTIONAL, INTENT(IN) :: ipType
!! interpolation type
!! used when baseInterpolation is Lagrange
INTEGER(I4B), OPTIONAL, INTENT(IN) :: basisType(:)
!! type of basis function used for
!! constructing the Lagrange polynomial
!! Used when baseInterpolation is Lagrange
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha(:)
!! alpha parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: beta(:)
!! beta parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda(:)
!! lambda parameter for Ultraspherical parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Ultraspherical
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: islocal
!! islocal denotes whether the order(:) is based on
!! local element or global element number.
!! local element means in order(ii) ii is the local
!! element number, global element means in order(ii) ii is the
!! global element number. Note that getting local element
!! number is difficult for user, so it is better to use
!! global element number.
END SUBROUTINE Initiate
END INTERFACE

Description

This subroutine configures a FEDOF object using inhomogeneous orders across mesh elements, allowing different polynomial orders for different parts of the domain.

Parameters

ParameterTypeIntentDescription
objCLASS(FEDOF_)INOUTThe Finite Element Degree of Freedom object to be initialized.
orderINTEGER(I4B)(:)INArray of polynomial orders for each cell element.
meshCLASS(AbstractMesh_), TARGETINThe computational mesh defining the spatial discretization.
baseContinuityCHARACTER(*)INSpecifies the continuity of basis functions (regularity).
baseInterpolationCHARACTER(*)INType of basis function used for interpolation.

Optional Parameters

ParameterTypeIntentDescription
ipTypeINTEGER(I4B)INInterpolation point type. Used when baseInterpolation is set to Lagrange.
basisTypeINTEGER(I4B)(:)INArray specifying types of basis functions used for constructing Lagrange polynomials. Used when baseInterpolation is Lagrange.
alphaREAL(DFP)(:)INAlpha parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
betaREAL(DFP)(:)INBeta parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
lambdaREAL(DFP)(:)INLambda parameters for Ultraspherical polynomials. Required when baseInterpolation is Lagrange and basisType is Ultraspherical.
islocalLOGICAL(LGT)INSpecifies whether order array references local element numbers (TRUE) or global element numbers (FALSE, default).

Usage Example

! Create an array with different orders for different elements
INTEGER(I4B) :: elementOrders(mesh%getTotalElements())
elementOrders = [1, 2, 2, 3, 2, 1] ! Example orders

CALL myFEDOF%Initiate(order=elementOrders, mesh=myMesh, &
baseContinuity="H1", baseInterpolation="Lagrange", &
ipType=GAUSS_LOBATTO, islocal=.FALSE.)

Notes

  • The length of the order array must match the number of elements in the mesh.
  • When islocal is not provided or is FALSE, the indices in order correspond to global element numbers.
  • When islocal is TRUE, the indices in order correspond to local element numbers, which may be different from global numbering.
  • Using global element numbering is generally easier for users to work with.

Interface 3

  • This method is used to initiate FEDOF by using ParameterList.
INTERFACE
MODULE SUBROUTINE Initiate(obj, param, mesh)
CLASS(FEDOF_), INTENT(INOUT) :: obj
TYPE(ParameterList_), INTENT(IN) :: param
CLASS(AbstractMesh_), TARGET, INTENT(IN) :: mesh
END SUBROUTINE Initiate
END INTERFACE
WIP

This interface is still under development and may not be fully functional yet.

Interface 4

  • This routine is similar to the Interface 2, but the order of the element is defined for global element numbers.
  • This method is more useful for the user who have no idea about the local element number.
  • order is a two-dimensional array.
    • The number of rows in order is equal to 2
    • The first row contains the global element number
    • The second row contains the order.
note

This routine will make order0(:) from order(:,:) and call Initiate2 method internally.

INTERFACE
MODULE SUBROUTINE obj_Initiate4(obj, order, mesh, baseContinuity, &
baseInterpolation, ipType, basisType, alpha, beta, lambda)
CLASS(FEDOF_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: order(:, :)
!! the number of columns in order is equal to total number of elements
!! the number of rows in order is equal to 2
!! the first row contains the global element number
!! the second rows contains the order of that element
CLASS(AbstractMesh_), TARGET, INTENT(IN) :: mesh
!! mesh
CHARACTER(*), INTENT(IN) :: baseContinuity
!! continuity of basis function
CHARACTER(*), INTENT(IN) :: baseInterpolation
!! interpolation of basis
INTEGER(I4B), OPTIONAL, INTENT(IN) :: ipType
!! interpolation type
!! used when baseInterpolation is Lagrange
INTEGER(I4B), OPTIONAL, INTENT(IN) :: basisType(:)
!! type of basis function used for
!! constructing the Lagrange polynomial
!! Used when baseInterpolation is Lagrange
REAL(DFP), OPTIONAL, INTENT(IN) :: alpha(:)
!! alpha parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: beta(:)
!! beta parameter for jacobian parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Jacobi
REAL(DFP), OPTIONAL, INTENT(IN) :: lambda(:)
!! lambda parameter for Ultraspherical parameter
!! used when baseInterpolation is Lagrange
!! used when basistype is Ultraspherical
END SUBROUTINE obj_Initiate4
END INTERFACE

Parameters

ParameterTypeIntentDescription
objCLASS(FEDOF_)INOUTThe FEDOF object to be initialized.
orderINTEGER(I4B)(:,:)IN2×N array where N is the number of elements to be configured.
meshCLASS(AbstractMesh_), TARGETINThe computational mesh defining the spatial discretization.
baseContinuityCHARACTER(*)INSpecifies the continuity of basis functions (regularity).
baseInterpolationCHARACTER(*)INType of basis function used for interpolation.

Optional Parameters

ParameterTypeIntentDescription
ipTypeINTEGER(I4B)INInterpolation point type. Used when baseInterpolation is set to Lagrange.
basisTypeINTEGER(I4B)(:)INArray specifying types of basis functions used for constructing Lagrange polynomials. Used when baseInterpolation is Lagrange.
alphaREAL(DFP)(:)INAlpha parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
betaREAL(DFP)(:)INBeta parameters for Jacobi polynomials. Required when baseInterpolation is Lagrange and basisType is Jacobi.
lambdaREAL(DFP)(:)INLambda parameters for Ultraspherical polynomials. Required when baseInterpolation is Lagrange and basisType is Ultraspherical.

Usage Example

! Create a 2×3 array to specify orders for three specific elements
INTEGER(I4B) :: elementOrders(2, 3)
! First row: global element numbers
elementOrders(1, :) = [1, 5, 10]
! Second row: corresponding polynomial orders
elementOrders(2, :) = [2, 3, 1]

CALL myFEDOF%Initiate(order=elementOrders, mesh=myMesh, &
baseContinuity="H1", baseInterpolation="Lagrange", &
ipType=GAUSS_LOBATTO)

Notes

  • This interface is more user-friendly as it allows specifying orders only for elements of interest using their global numbers.
  • Internally, this method constructs a complete order array and calls Interface 2.
  • The number of columns in the order array can be less than the total number of elements in the mesh - only specified elements will receive custom orders.

Example (H1, Hierarchical, Uniform Order)

This example shows how to initiate FEDOF with H1 and Hierarchy interpolation.

!> author: Vikas Sharma, Ph. D.
! date: 2025-06-06
! summary: Test the intitiate method for H1, Heirarchical basis,
! for different orders.

PROGRAM main
USE FEDOF_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE HDF5File_Class
USE Display_Method
USE GlobalData
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION

IMPLICIT NONE

TYPE(FEDOF_) :: obj
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
CHARACTER(*), PARAMETER :: filename = &
"../../FEMesh/examples/meshdata/small_tri3_mesh.h5"
TYPE(HDF5File_) :: meshfile
INTEGER(I4B) :: found, want

CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL meshfile%Initiate(filename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(meshfile, '')

meshptr => dom%GetMeshPointer()

CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Heirarchical", &
order=1, mesh=meshptr)
!CALL fedof%Display("FEDOF:")
found = obj%GetTotalDOF()
want = meshptr%GetTotalNodes()
CALL OK(found == want, "Total DOF (order=1): ")

CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Heirarchical", &
order=2, mesh=meshptr)
found = obj%GetTotalDOF()
want = meshptr%GetTotalNodes() + meshptr%GetTotalFaces()
CALL OK(found == want, "Total DOF (order=2): ")

CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Heirarchical", &
order=3, mesh=meshptr)
found = obj%GetTotalDOF()
want = meshptr%GetTotalNodes() + 2*meshptr%GetTotalFaces() + meshptr%GetTotalCells()
CALL OK(found == want, "Total DOF (order=3): ")

CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Heirarchical", &
order=4, mesh=meshptr)
found = obj%GetTotalDOF()
want = meshptr%GetTotalNodes() + 3*meshptr%GetTotalFaces() + 3*meshptr%GetTotalCells()
CALL OK(found == want, "Total DOF (order=4): ")

!CALL dom%Display("domain:")
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()

END PROGRAM main

Example (H1, Hierarchical, Inhomogeneous Order)

This example shows how to initiate FEDOF with H1 and Hierarchy interpolation. Interface 2 of initiate method is used wherein order is a vector of integers representing the order of basis functions for each cell.

!> author: Vikas Sharma, Ph. D.
! date: 2025-06-01
! summary: Initiate fedof for H1 and Heirarchical bases, order is a vector.

PROGRAM main
USE FEDOF_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE HDF5File_Class
USE Display_Method
USE GlobalData
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ReallocateUtility

IMPLICIT NONE

TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
CHARACTER(*), PARAMETER :: &
filename = "../../FEMesh/examples/meshdata/small_tri3_mesh.h5", &
baseInterpolation = "Hierarchical", &
baseContinuity = "H1"

TYPE(HDF5File_) :: meshfile
INTEGER(I4B) :: found, want
INTEGER(I4B), ALLOCATABLE :: cellOrder(:)

CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL meshfile%Initiate(filename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(meshfile, '')

meshptr => dom%GetMeshPointer()
! CALL meshptr%DisplayMeshInfo("Mesh Info:")

CALL Reallocate(cellOrder, meshptr%GetTotalCells())
cellOrder = 1

CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=cellOrder, mesh=meshptr, islocal=.TRUE.)
! CALL fedof%Display("FEDOF:")
found = fedof%GetTotalDOF()
want = meshptr%GetTotalNodes()
CALL IS(found, want, "Total DOF (order=1): ")

cellOrder = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=cellOrder, mesh=meshptr, islocal=.TRUE.)

found = fedof%GetTotalDOF()
want = meshptr%GetTotalNodes() + meshptr%GetTotalFaces()
CALL IS(found, want, "Total DOF (order=2): ")

cellOrder = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=cellOrder, mesh=meshptr, islocal=.TRUE.)

found = fedof%GetTotalDOF()
want = meshptr%GetTotalNodes() + 2*meshptr%GetTotalFaces() + meshptr%GetTotalCells()
CALL IS(found, want, "Total DOF (order=3): ")

cellOrder = 4
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=cellOrder, mesh=meshptr, islocal=.TRUE.)
found = fedof%GetTotalDOF()
want = meshptr%GetTotalNodes() + 3*meshptr%GetTotalFaces() + 3*meshptr%GetTotalCells()
CALL IS(found, want, "Total DOF (order=4): ")

!CALL dom%Display("domain:")
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()

END PROGRAM main

Example (H1, Hierarchical, Inhomogeneous Order with Local Element Number)

!> author: Vikas Sharma, Ph. D.
! date: 2025-06-06
! summary: Initiate fedof with H1 and Heirarchical bases, order is a vector.

PROGRAM main
USE FEDOF_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE HDF5File_Class
USE Display_Method
USE GlobalData
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ReallocateUtility

IMPLICIT NONE

TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
CHARACTER(*), PARAMETER :: &
filename = "../../FEMesh/examples/meshdata/small_tri3_mesh.h5", &
baseInterpolation = "Hierarchical", &
baseContinuity = "H1"
TYPE(HDF5File_) :: meshfile
LOGICAL(LGT) :: isok
INTEGER(I4B) :: found, want, order, ii, iel
INTEGER(I4B), ALLOCATABLE :: cellOrder(:), aintvec(:)

CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL meshfile%Initiate(filename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(meshfile, '')

meshptr => dom%GetMeshPointer()
CALL Reallocate(cellOrder, meshptr%GetTotalCells())

aintvec = [14, 15, 17, 24, 18]
order = 1
DO iel = 1, SIZE(aintvec)
isok = meshptr%IsElementPresent(aintvec(iel))
IF (.NOT. isok) CYCLE

ii = meshptr%GetLocalElemNumber(aintvec(iel))
cellOrder(ii) = order
END DO

aintvec = [22, 23, 25, 21]
order = 2
DO iel = 1, SIZE(aintvec)
isok = meshptr%IsElementPresent(aintvec(iel))
IF (.NOT. isok) CYCLE

ii = meshptr%GetLocalElemNumber(aintvec(iel))
cellOrder(ii) = order
END DO

aintvec = [19, 20, 26, 16, 13]
order = 3
DO iel = 1, SIZE(aintvec)
isok = meshptr%IsElementPresent(aintvec(iel))
IF (.NOT. isok) CYCLE

ii = meshptr%GetLocalElemNumber(aintvec(iel))
cellOrder(ii) = order
END DO

CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=cellOrder, mesh=meshptr, islocal=.TRUE.)
found = fedof%GetTotalDOF()
want = 39
isok = found == want
CALL OK(isok, "Total DOF ")
IF (.NOT. isok) CALL Display([found, want], "found, want: ")

CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()

END PROGRAM main

Example (H1, Lagrange, Uniform Order)

!> author: Vikas Sharma, Ph. D.
! date: 2024-05-24
! summary: Lagrange polynomial is tested in this example

PROGRAM main
USE FEDOF_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE HDF5File_Class
USE Display_Method
USE GlobalData
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
use BaseType, only: poly=>TypePolynomialOpt

IMPLICIT NONE

TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
CHARACTER(*), PARAMETER :: filename = &
"../../FEMesh/examples/meshdata/small_tri3_mesh.h5"

TYPE(HDF5File_) :: meshfile
INTEGER(I4B) :: found, want
INTEGER(I4B), PARAMETER :: order = 1, ipType = poly%monomial
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Lagrange"

CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL meshfile%Initiate(filename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(meshfile, '')

meshptr => dom%GetMeshPointer()

CALL fedof%Initiate(baseContinuity=baseContinuity, ipType=ipType, &
baseInterpolation=baseInterpolation, order=order, mesh=meshptr)
!CALL fedof%Display("FEDOF:")
found = fedof%GetTotalDOF()
want = meshptr%GetTotalNodes()
CALL IS(found, want, "Total DOF (order=1): ")

!CALL dom%Display("domain:")
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()

END PROGRAM main

Example (H1, Lagrange, Inhomogeneous Order)