Skip to main content

Copy

The Copy method copies the contents of one FEDOF_ object to another. This method is also used by the assignment operator (=) for FEDOF_ objects.

Interface

INTERFACE
MODULE SUBROUTINE Copy(obj, obj2)
CLASS(FEDOF_), INTENT(INOUT) :: obj
CLASS(FEDOF_), INTENT(IN) :: obj2
END SUBROUTINE Copy
END INTERFACE

Parameters

  • obj - The destination FEDOF_ object that will receive the copied data
  • obj2 - The source FEDOF_ object whose data will be copied

Implementation Details

The method performs a shallow copy of most attributes from obj2 to obj:

  1. Copies basic properties:

    • isLagrange - Flag indicating if the base interpolation is Lagrange
    • tdof - Total number of degrees of freedom
    • tNodes, tEdges, tFaces, tCells - Total count of different mesh entities
    • baseContinuity - Continuity or conformity of basis (e.g., "H1")
    • baseInterpolation - Type of basis functions (e.g., "LAGR", "HIER")
    • maxCellOrder, maxFaceOrder, maxEdgeOrder - Maximum order values
  2. Sets the mesh pointer to point to the same mesh as the source object:

obj%mesh => obj2%mesh
  1. Copies array data when allocated in the source object:

    • cellOrder - Order of each cell
    • faceOrder - Order of each face
    • edgeOrder - Order of each edge
    • edgeIA, faceIA, cellIA - Sparsity patterns for different entities
  2. Associates finite element pointers:

DO ii = 1, SIZE(obj2%fe)
isok = ASSOCIATED(obj2%fe(ii)%ptr)
IF (isok) THEN
obj%fe(ii)%ptr => obj2%fe(ii)%ptr
END IF
END DO

Notes

  • This is a shallow copy - the mesh and finite element pointers in the destination object point to the same objects as the source.
  • Array data is copied only if it exists in the source object.
  • The method doesn't allocate new memory for arrays; it assumes the destination arrays either don't exist or are already properly allocated.
  • This method is particularly useful when you need to create a duplicate of a FEDOF_ object without creating new underlying mesh or finite element instances.

Example

!> author: Vikas Sharma, Ph. D.
! date: 2025-06-06
! summary: This method tests the COPY method for FEDOF class.

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

IMPLICIT NONE

TYPE(FEDOF_) :: obj, obj2
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: meshptr => NULL()
CHARACTER(*), PARAMETER :: &
filename = "../../FEMesh/examples/meshdata/small_tri3_mesh.h5", &
baseInterpolation = "Heirarchical", &
baseContinuity = "H1"

TYPE(HDF5File_) :: meshfile
INTEGER(I4B) :: 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=baseContinuity, &
baseInterpolation=baseInterpolation, order=1, mesh=meshptr)
!CALL fedof%Display("FEDOF:")
CALL obj2%copy(obj)

found = obj2%GetTotalDOF()
want = meshptr%GetTotalNodes()
CALL OK(found == want, "Total DOF (order=1): ")

CALL obj%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=2, mesh=meshptr)
obj2 = obj
found = obj2%GetTotalDOF()
want = meshptr%GetTotalNodes() + meshptr%GetTotalFaces()
CALL OK(found == want, "Total DOF (order=2): ")

CALL obj%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=3, mesh=meshptr)
CALL obj2%Copy(obj)

found = obj2%GetTotalDOF()
want = meshptr%GetTotalNodes() + 2*meshptr%GetTotalFaces() + meshptr%GetTotalCells()
CALL OK(found == want, "Total DOF (order=3): ")

CALL obj%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=4, mesh=meshptr)
obj2 = obj

found = obj2%GetTotalDOF()
want = meshptr%GetTotalNodes() + 3*meshptr%GetTotalFaces() + 3*meshptr%GetTotalCells()
CALL OK(found == want, "Total DOF (order=4): ")

!CALL dom%Display("domain:")
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()

END PROGRAM main