GetCaseName
The GetCaseName
method returns a string that identifies the combination of continuity and interpolation type used in the finite element degree of freedom object. This case name is used internally to determine which algorithms to apply when processing the finite element data.
Interface
INTERFACE
MODULE FUNCTION GetCaseName(obj) RESULT(ans)
CLASS(FEDOF_), INTENT(IN) :: obj
CHARACTER(6) :: ans
END FUNCTION GetCaseName
END INTERFACE
Parameters
obj
- TheFEDOF_
object from which to obtain the case name
Return Value
ans
- A 6-character string combining the basis continuity and the interpolation type
Implementation Details
The method is straightforward, concatenating the baseContinuity
and baseInterpolation
properties of the FEDOF_
object:
ans = obj%baseContinuity//obj%baseInterpolation
Possible Return Values
The return value is a combination of the continuity type and interpolation type:
-
For continuity (
obj%baseContinuity
):"H1"
- Standard conforming elements"HC[url]"
- H(curl) conforming elements"HD[iv]"
- H(div) conforming elements"DG"
- Discontinuous Galerkin elements
-
For interpolation (
obj%baseInterpolation
):"LAGR"
- Lagrange interpolation"HIER"
- Hierarchical interpolation"ORTHO"
- Orthogonal interpolation"HERM"
- Hermite interpolation"SERE"
- Serendipity interpolation
Common combinations include:
"H1LAGR"
- H¹ conforming elements with Lagrange interpolation"H1HIER"
- H¹ conforming elements with hierarchical interpolation"DGHIER"
- Discontinuous Galerkin with hierarchical interpolation
Usage Examples
The case name is primarily used in conditional branching to select the appropriate algorithms for a given finite element type:
CHARACTER(6) :: casename
casename = fedof%GetCaseName()
SELECT CASE (casename)
CASE ('H1LAGR')
! Process H¹ conforming elements with Lagrange interpolation
CASE ('H1HIER')
! Process H¹ conforming elements with hierarchical interpolation
CASE ('DGHIER')
! Process Discontinuous Galerkin with hierarchical interpolation
CASE DEFAULT
! Handle unexpected case
END SELECT
This method is used extensively within the FEDOF_
class to determine the appropriate algorithms for shape function calculation, connectivity generation, and other operations that depend on the specific type of finite element formulation.
Example
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Get case name
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: TypeQuadratureOpt
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
CHARACTER(6) :: 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)
found = obj%GetCaseName()
want = "H1HEIR"
CALL OK(found == want, "GetCaseName H1HEIR: ")
CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Hierarchical", &
order=1, mesh=meshptr)
found = obj%GetCaseName()
want = "H1HIER"
CALL OK(found == want, "GetCaseName H1HIER: ")
CALL obj%Initiate(baseContinuity="H1", baseInterpolation="Lagrange", &
order=1, mesh=meshptr, ipType=TypeQuadratureOpt%equidistance)
found = obj%GetCaseName()
want = "H1LAGR"
CALL OK(found == want, "GetCaseName H1LAGR: ")
!CALL dom%Display("domain:")
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
END PROGRAM main