GetFaceDOF
The GetFaceDOF
generic method in the FEDOF_Class
is used to retrieve the degrees of freedom (DOF) associated with a face in a finite element mesh. This method has two specific implementations with different parameter sets.
Method Description
The GetFaceDOF
generic method retrieves face degrees of freedom through two different interfaces:
GetFaceDOF1
- Gets DOFs directly using a global face numberGetFaceDOF2
- Gets DOFs using a global element number and local face number within that element
Interface 1
SUBROUTINE obj_GetFaceDOF1(obj, globalFace, ans, tsize, islocal)
Parameters
obj
- Input,CLASS(FEDOF_)
, FEDOF object instanceglobalFace
- Input,INTEGER(I4B)
, global face numberans
- Output,INTEGER(I4B)(:)
, array to store face degrees of freedomtsize
- Output,INTEGER(I4B)
, total size of data written toans
islocal
- Optional Input,LOGICAL(LGT)
, if true,globalFace
is treated as a local face number
Implementation Details
The implementation in GetMethods
submodule uses the internal face index array structure:
INTEGER(I4B) :: ii
tsize = 0
DO ii = obj%faceIA(globalface), obj%faceIA(globalface + 1) - 1
tsize = tsize + 1
ans(tsize) = ii
END DO
The method retrieves DOFs from the face sparsity data structure using the compressed sparse row format stored in faceIA
.
Interface 2
SUBROUTINE obj_GetFaceDOF2(obj, globalElement, localFaceNumber, ans, tsize, islocal)
Parameters
obj
- Input,CLASS(FEDOF_)
, DOF objectglobalElement
- Input,INTEGER(I4B)
, global or local element numberlocalFaceNumber
- Input,INTEGER(I4B)
, local face number in global elementans
- Output,INTEGER(I4B)(:)
, array to store face degrees of freedomtsize
- Output,INTEGER(I4B)
, total size of data written toans
islocal
- Optional Input,LOGICAL(LGT)
, if true,globalElement
is treated as a local element number
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)
CALL obj%GetFaceDOF(globalFace=globalFace, ans=ans, tsize=tsize, &
islocal=islocal)
Usage Example
! Example to get DOFs from a face
INTEGER(I4B) :: faceDOFs(100), totalDOFs
TYPE(FEDOF_) :: myDOF
! Method 1: Using global face number
CALL myDOF%GetFaceDOF(globalFace=5, ans=faceDOFs, tsize=totalDOFs)
! Method 2: Using element and local face number
CALL myDOF%GetFaceDOF(globalElement=10, localFaceNumber=2, ans=faceDOFs, tsize=totalDOFs)
Related Methods
GetTotalFaceDOF
- Returns the total number of DOFs on a faceGetEdgeDOF
- Similar function but for element edgesGetCellDOF
- Gets DOFs for volumetric elementsGetVertexDOF
- Gets DOFs for element vertices
Need any additional details about these methods?
Example
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Testing GetFaceDOF 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//" GetFaceDOF 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(10), want(10)
INTEGER(I4B), ALLOCATABLE :: con(:)
order = 1
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
con = fedof%GetConnectivity(globalElement=1, islocal=.TRUE., opt="A")
CALL fedof%GetFaceDOF(globalFace=1, islocal=.TRUE., &
ans=found, tsize=tsize)
isok = tsize .EQ. 0
CALL OK(isok, testname//" (order= "//ToString(order)//"): ")
CALL fedof%GetFaceDOF(globalElement=1, islocal=.TRUE., &
ans=found, tsize=tsize, localFaceNumber=1)
isok = tsize .EQ. 0
CALL OK(isok, testname//" interface 2 (order= "//ToString(order)//"): ")
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
INTEGER(I4B) :: tsize, found(10), want(10)
order = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL fedof%GetFaceDOF(globalFace=1, islocal=.TRUE., &
ans=found, tsize=tsize)
CALL IS(tsize, 1, testname//" (order= "//ToString(order)//"): ")
CALL IS(found(1), 13, testname//" (order= "//ToString(order)//"): ")
CALL fedof%GetFaceDOF(globalElement=1, islocal=.TRUE., &
ans=found, tsize=tsize, localFaceNumber=1)
CALL IS(tsize, 1, testname//" (order= "//ToString(order)//"): ")
CALL IS(found(1), 13, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test3
INTEGER(I4B) :: tsize, found(10), want(10)
order = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL fedof%GetFaceDOF(globalFace=1, islocal=.TRUE., &
ans=found, tsize=tsize)
CALL IS(tsize, 2, testname//" (order= "//ToString(order)//"): ")
CALL IS(found(1), 13, testname//" (order= "//ToString(order)//"): ")
CALL fedof%GetFaceDOF(globalElement=1, islocal=.TRUE., &
ans=found, tsize=tsize, localFaceNumber=1)
CALL IS(tsize, 2, testname//" (order= "//ToString(order)//"): ")
CALL IS(found(1), 13, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test3
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main