GetTotalCellDOF
The GetTotalCellDOF
method in the FEDOF_Class
provides the total number of degrees of freedom (DOF) associated with a cell in a finite element mesh.
Interface
From the FEDOF_Class.F90
file, the GetTotalCellDOF
method has the following interface:
PROCEDURE, PUBLIC, PASS(obj) :: GetTotalCellDOF => obj_GetTotalCellDOF
!! Get total cell degrees of freedom
The specific interface is defined as:
MODULE FUNCTION obj_GetTotalCellDOF(obj, globalCell, islocal) RESULT(ans)
CLASS(FEDOF_), INTENT(IN) :: obj
INTEGER(I4B), INTENT(IN) :: globalCell
LOGICAL(LGT), INTENT(IN), OPTIONAL :: islocal
INTEGER(I4B) :: ans
END FUNCTION obj_GetTotalCellDOF
Description
The GetTotalCellDOF
method calculates the total number of degrees of freedom associated with a specific cell (volumetric element) in a finite element mesh.
Parameters
obj
- Input,CLASS(FEDOF_)
, FEDOF object instanceglobalCell
- Input,INTEGER(I4B)
, global cell number or local cell number (depends on islocal)islocal
- Optional Input,LOGICAL(LGT)
, if true,globalCell
is treated as a local cell numberans
- Output,INTEGER(I4B)
, total number of DOFs on the specified cell
Implementation Details
The implementation uses the internal cell index array structure to calculate the total DOF count:
INTEGER(I4B) :: jj
jj = obj%mesh%GetLocalElemNumber(globalElement=globalCell, islocal=islocal)
ans = obj%cellIA(jj + 1) - obj%cellIA(jj)
This efficiently computes the number of DOFs by finding the difference between consecutive indices in the compressed sparse row format that stores the cell degrees of freedom.
Usage Example
! Example to get total DOFs on a cell
INTEGER(I4B) :: totalDOFs
TYPE(FEDOF_) :: myDOF
! Get total DOFs on cell #5
totalDOFs = myDOF%GetTotalCellDOF(globalCell=5)
! Get total DOFs using local cell numbering
totalDOFs = myDOF%GetTotalCellDOF(globalCell=2, islocal=.TRUE.)
Related Methods
GetCellDOF
- Retrieves the actual DOF indices for a cellGetTotalEdgeDOF
- Gets total DOF count for element edgesGetTotalFaceDOF
- Gets total DOF count for element facesGetTotalVertexDOF
- Gets total DOF count for verticesGetTotalDOF
- Gets total DOF count for an entire element or the entire mesh
These methods together provide a complete interface for accessing and counting degrees of freedom in different components of the finite element mesh hierarchy (vertices, edges, faces, cells).
Example 1
!> author: Vikas Sharma, Ph. D.
! date: 2025-06-07
! summary: Testing GetTotalCellDOF 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//" GetTotalCellDOF 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 test4
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
INTEGER(I4B) :: tsize, want, found
order = 1
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalCellDOF(globalCell=1, islocal=.TRUE.)
want = 0
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
INTEGER(I4B) :: tsize, want, found
order = 2
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalCellDOF(globalCell=1, islocal=.TRUE.)
want = 0
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test2
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test3
INTEGER(I4B) :: tsize, want, found
order = 3
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalCellDOF(globalCell=1, islocal=.TRUE.)
want = 1
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test3
!----------------------------------------------------------------------------
! test3
!----------------------------------------------------------------------------
SUBROUTINE test4
INTEGER(I4B) :: tsize, want, found
order = 4
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, &
order=order, mesh=meshptr)
found = fedof%GetTotalCellDOF(globalCell=1, islocal=.TRUE.)
want = 3
CALL IS(found, want, testname//" (order= "//ToString(order)//"): ")
END SUBROUTINE test4
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main