GetTotalFaceDOF
The GetTotalFaceDOF
generic method in the FEDOF_Class
provides the total number of degrees of freedom (DOF) associated with a face in a finite element mesh.
Interfaces
From the FEDOF_Class.F90
file, the GetTotalFaceDOF
generic method has the following interfaces:
PROCEDURE, PASS(obj) :: GetTotalFaceDOF1 => GetTotalFaceDOF1
!! Get total face dof
PROCEDURE, PASS(obj) :: GetTotalFaceDOF2 => GetTotalFaceDOF2
!! Get total face dof from global element and local face number
GENERIC, PUBLIC :: GetTotalFaceDOF => GetTotalFaceDOF1, GetTotalFaceDOF2
The specific interfaces are defined as:
MODULE FUNCTION GetTotalFaceDOF1(obj, globalFace, islocal) RESULT(ans)
CLASS(FEDOF_), INTENT(IN) :: obj
INTEGER(I4B), INTENT(IN) :: globalFace
LOGICAL(LGT), INTENT(IN), OPTIONAL :: islocal
INTEGER(I4B) :: ans
END FUNCTION GetTotalFaceDOF1
MODULE FUNCTION GetTotalFaceDOF2(obj, globalElement, localFaceNumber, &
islocal) RESULT(ans)
CLASS(FEDOF_), INTENT(IN) :: obj
!! DOF object
INTEGER(I4B), INTENT(IN) :: globalElement
!! global or local element number
INTEGER(I4B), INTENT(IN) :: localFaceNumber
!! local face number in globall element
LOGICAL(LGT), INTENT(IN), OPTIONAL :: islocal
!! if true then globalElement is local element
INTEGER(I4B) :: ans
!! Total number of degree of freedom on face
END FUNCTION GetTotalFaceDOF2
Description
The GetTotalFaceDOF
generic method provides two implementations to retrieve the total face DOF count:
GetTotalFaceDOF1
- Gets total DOF count directly using a global face numberGetTotalFaceDOF2
- Gets total DOF count using a global element number and local face number within that element
Interface 1: GetTotalFaceDOF1
Parameters
obj
- Input,CLASS(FEDOF_)
, FEDOF object instanceglobalFace
- Input,INTEGER(I4B)
, global face numberislocal
- Optional Input,LOGICAL(LGT)
, if true,globalFace
is treated as a local face numberans
- Output,INTEGER(I4B)
, total number of DOFs on the specified face
Implementation Details
The implementation uses the internal face index array structure to calculate the total DOF count:
ans = obj%faceIA(globalface + 1) - obj%faceIA(globalface)
This efficiently computes the number of DOFs by finding the difference between consecutive indices in the compressed sparse row format.
Interface 2: GetTotalFaceDOF2
Parameters
obj
- Input,CLASS(FEDOF_)
, DOF objectglobalElement
- Input,INTEGER(I4B)
, global or local element numberlocalFaceNumber
- Input,INTEGER(I4B)
, local face number in global elementislocal
- Optional Input,LOGICAL(LGT)
, if true,globalElement
is treated as a local elementans
- Output,INTEGER(I4B)
, total number of degrees of freedom on the specified face
Implementation Details
This implementation first converts the element and local face information to a global face number, then calls the first implementation:
INTEGER(I4B) :: globalFace
globalFace = obj%mesh%GetGlobalFaceNumber(globalElement=globalElement, &
islocal=islocal, localFaceNumber=localFaceNumber)
ans = obj%GetTotalFaceDOF(globalFace=globalFace, islocal=islocal)
Usage Example
! Example to get total DOFs on a face
INTEGER(I4B) :: totalDOFs
TYPE(FEDOF_) :: myDOF
! Method 1: Using global face number
totalDOFs = myDOF%GetTotalFaceDOF(globalFace=3)
! Method 2: Using element and local face number
totalDOFs = myDOF%GetTotalFaceDOF(globalElement=7, localFaceNumber=2)
Related Methods
GetFaceDOF
- Retrieves the actual DOF indices for a faceGetTotalEdgeDOF
- Similar function but for element edgesGetTotalDOF
- Gets total DOF count for an entire element or the entire meshGetTotalVertexDOF
- Gets total DOF count for vertices
These methods are essential for traversing and accessing degrees of freedom in different components of the mesh hierarchy (vertices, edges, faces, cells).
Example 1
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Testing GetTotalFaceDOF method of FEDOF class
! H1 Heirarchical Second Order Triangular Mesh
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 AppendUtility
USE ArangeUtility
USE ReallocateUtility
IMPLICIT NONE
CHARACTER(*), PARAMETER :: &
filename = "../../FEMesh/examples/meshdata/small_tri3_mesh.h5", &
baseContinuity = "H1", &
baseInterpolation = "Heirarchical", &
testname = baseContinuity//" "//baseInterpolation//" GetTotalFaceDOF test"
TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
TYPE(HDF5File_) :: meshfile
INTEGER(I4B) :: order, totalVertexNodes, totalFaces
LOGICAL(LGT) :: isok
CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL meshfile%Initiate(filename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(meshfile, '')
meshptr => dom%GetMeshPointer()
totalVertexNodes = meshptr%GetTotalVertexNodes()
totalFaces = meshptr%GetTotalFaces()
CALL test1
CALL test2
CALL test3
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
INTEGER(I4B) :: tsize, found, want
INTEGER(I4B), ALLOCATABLE :: con(:)
order = 1
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalFaceDOF(globalFace=1, islocal=.TRUE.)
want = 0
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
found = fedof%GetTotalFaceDOF(globalElement=1, islocal=.TRUE., &
localFaceNumber=1)
CALL IS(found, want, testname//" interface 2 (order= "//ToString(order)//"): ")
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
INTEGER(I4B) :: tsize, found, want
order = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalFaceDOF(globalFace=1, islocal=.TRUE.)
want = 1
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
found = fedof%GetTotalFaceDOF(globalElement=1, islocal=.TRUE., &
localFaceNumber=1)
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test3
INTEGER(I4B) :: tsize, found, want
order = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalFaceDOF(globalFace=1, islocal=.TRUE.)
want = 2
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
found = fedof%GetTotalFaceDOF(globalElement=1, islocal=.TRUE., &
localFaceNumber=1)
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test3
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main