GetVertexDOF
The GetVertexDOF
method retrieves the degree of freedom (DOF) number associated with a specified vertex (node) in the mesh. In finite element analysis, each vertex in the mesh is typically assigned one or more degrees of freedom, which represent the unknown values to be solved for.
Interface
INTERFACE
MODULE SUBROUTINE GetVertexDOF(obj, globalNode, ans, tsize, islocal)
CLASS(FEDOF_), INTENT(IN) :: obj
INTEGER(I4B), INTENT(IN) :: globalNode
INTEGER(I4B), INTENT(INOUT) :: ans(:)
INTEGER(I4B), INTENT(OUT) :: tsize
LOGICAL(LGT), INTENT(IN), OPTIONAL :: islocal
END SUBROUTINE GetVertexDOF
END INTERFACE
Parameters
obj
- TheFEDOF_
object containing degree of freedom informationglobalNode
- The global or local node number, depending on the value ofislocal
ans
- An array that will store the degree of freedom number(s) associated with the specified nodetsize
- The number of degrees of freedom written to theans
arrayislocal
- Optional logical flag. If present and.TRUE.
,globalNode
is interpreted as a local node number; otherwise, it's a global node number
Implementation Details
The implementation is straightforward for vertex DOFs:
tsize = 1
ans(1) = obj%mesh%GetLocalNodeNumber(globalNode, islocal=islocal)
For vertex DOFs in most finite element formulations:
- The method sets
tsize
to 1, indicating that one degree of freedom is associated with each vertex - It obtains the local node number using the mesh's
GetLocalNodeNumber
method, which converts global node numbers to local ones if necessary - The local node number itself is used as the degree of freedom number
Notes
- This method assumes a simple one-to-one mapping between mesh vertices and degrees of freedom, which is typical for scalar problems (e.g., heat conduction, potential flow)
- For vector problems (e.g., elasticity) or higher-order problems, a more complex mapping might be needed, potentially requiring multiple DOFs per vertex
- The implementation reflects that in typical finite element numbering schemes, vertex DOFs are numbered first, followed by edge, face, and cell DOFs
- In the current implementation, the DOF number for a vertex is identical to its local node number in the mesh
Example Usage
This method is typically used when assembling the global system of equations or when mapping local element matrices to the global system:
INTEGER(I4B) :: nodeDOF(1), tsize, nodeNum
! Get the DOF number for node 5
nodeNum = 5
CALL fedof%GetVertexDOF(globalNode=nodeNum, ans=nodeDOF, tsize=tsize)
! nodeDOF(1) now contains the DOF number for node 5
! tsize will be 1
This method is often called internally by other methods in the FEDOF_
class, particularly when constructing element connectivity arrays that map local element degrees of freedom to global system degrees of freedom.
Example 1
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Testing GetVertexDOF 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
TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
TYPE(HDF5File_) :: meshfile
INTEGER(I4B), ALLOCATABLE :: found(:), want(:)
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
order = 1
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF (order= "//ToString(order)//"): ")
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
INTEGER(I4B) :: tsize
order = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF (order= "//ToString(order)//"): ")
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test3
INTEGER(I4B) :: tsize
order = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF (order= "//ToString(order)//"): ")
END SUBROUTINE test3
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main
Example 2
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Testing GetVertexDOF 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_tri6_mesh.h5", &
baseContinuity = "H1", &
baseInterpolation = "Heirarchical", &
testname = baseContinuity//" "//baseInterpolation
TYPE(FEDOF_) :: fedof
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
TYPE(HDF5File_) :: meshfile
INTEGER(I4B), ALLOCATABLE :: found(:), want(:)
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
order = 1
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF "//ToString(order)//"): ")
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
INTEGER(I4B) :: tsize
order = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF (order= "//ToString(order)//"): ")
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test3
INTEGER(I4B) :: tsize
order = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
CALL Reallocate(found, 2, want, 2)
CALL fedof%GetVertexDOF(globalNode=1, islocal=.TRUE., &
ans=found, tsize=tsize)
want(1) = 1
isok = ALL(found == want)
CALL OK(isok, testname//" GetVertexDOF (order= "//ToString(order)//"): ")
END SUBROUTINE test3
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main