Set
Set the entries in the ScalarField.
Calling example:
Set a single entry
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: globalNode
REAL(DFP):: VALUE
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Set all entries to a single scalar value
CALL Set(
CLASS(ScalarField_):: obj
REAL(DFP):: VALUE
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Set all entries of scalar field to a given vector.
CALL Set(
CLASS(ScalarField_):: obj
REAL(DFP):: VALUE(:)
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Set selected entries to a single value.
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: globalNode(:)
REAL(DFP):: VALUE
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Set multiple entries to different values.
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: globalNode(:)
REAL(DFP):: VALUE(:)
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Select multiple enties using triads.
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: istart
INTEGER(I4B):: iend
INTEGER(I4B):: stride
REAL(DFP):: VALUE
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Select multiple enties using triads.
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: istart
INTEGER(I4B):: iend
INTEGER(I4B):: stride
REAL(DFP):: VALUE(:)
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Copy obj2 in obj.
CALL Set(
CLASS(ScalarField_):: obj
CLASS(ScalarField_):: obj2
)
Select multiple values using FEVariable.
CALL Set(
CLASS(ScalarField_):: obj
INTEGER(I4B):: globalNode(:)
TYPE(FEVariable_):: VALUE
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Copy obj2 into obj like AXPY.
CALL Set(
CLASS(ScalarField_):: obj
CLASS(ScalarField_):: obj2
REAL(DFP):: scale
LOGICAL(LGT):: addContribution
)
Interface 1
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set1(obj, globalNode, VALUE, scale, &
& addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: globalNode
REAL(DFP), INTENT(IN) :: VALUE
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set1
END INTERFACE
- Set single entry.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set1, setting a single value
PROGRAM main
USE AbstractField_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE GlobalData
USE ScalarField_Class
USE FEDOF_Class
USE Display_Method
USE ApproxUtility
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
IMPLICIT NONE
CHARACTER(*), PARAMETER :: tomlFileName = "./Set1.toml", &
myName = "main", &
modName = "_ImportFromToml_test_1.F90"
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(FEDOF_) :: fedof, geofedof
TYPE(ScalarField_) :: u
CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL dom%ImportFromToml(fileName=tomlFileName, tomlName="domain")
CALL u%ImportFromToml(fedof=fedof, geofedof=geofedof, dom=dom, &
fileName=tomlFileName, tomlName="u")
mesh => dom%GetMeshPointer()
CALL test1
CALL test2
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
CHARACTER(*), PARAMETER :: testName = "test1()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans
LOGICAL(LGT) :: isok
tNodes = u%SIZE()
DO ii = 1, tNodes
areal = REAL(ii, DFP)
CALL u%Set(globalNode=ii, islocal=.TRUE., VALUE=areal)
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
isok = areal.APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
CHARACTER(*), PARAMETER :: testName = "test2()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans, scale
LOGICAL(LGT) :: isok
tNodes = u%SIZE()
scale = 2.0_DFP
DO ii = 1, tNodes
areal = REAL(ii, DFP)
CALL u%Set(globalNode=ii, islocal=.TRUE., VALUE=areal, scale=scale)
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
areal = areal * scale
isok = areal.APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test2
END PROGRAM main
Interface 2
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set2(obj, VALUE, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
REAL(DFP), INTENT(IN) :: VALUE
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set2
END INTERFACE
- Set all values of scalar field to a given scalar.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set2, setting all the values
PROGRAM main
USE AbstractField_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE GlobalData
USE ScalarField_Class
USE FEDOF_Class
USE Display_Method
USE ApproxUtility
USE Test_Method
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
IMPLICIT NONE
CHARACTER(*), PARAMETER :: tomlFileName = "./Set1.toml", &
myName = "main", &
modName = "_ImportFromToml_test_1.F90"
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(FEDOF_) :: fedof, geofedof
TYPE(ScalarField_) :: u
CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL dom%ImportFromToml(fileName=tomlFileName, tomlName="domain")
CALL u%ImportFromToml(fedof=fedof, geofedof=geofedof, dom=dom, &
fileName=tomlFileName, tomlName="u")
mesh => dom%GetMeshPointer()
CALL test1
CALL test2
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
CHARACTER(*), PARAMETER :: testName = "test1()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans
LOGICAL(LGT) :: isok
tNodes = u%SIZE()
areal = 5.0_DFP
CALL u%Set(VALUE=areal)
DO ii = 1, tNodes
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
isok = areal.APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
CHARACTER(*), PARAMETER :: testName = "test2()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans, scale
LOGICAL(LGT) :: isok
tNodes = u%SIZE()
scale = 2.0_DFP
areal = 5.0_DFP
CALL u%Set(VALUE=areal, scale=scale)
DO ii = 1, tNodes
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
isok = (areal * scale) .APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
CALL Display(" Expected: "//ToString(areal * scale)// &
" Found: "//ToString(ans))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test2
END PROGRAM main
Interface 3
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set3(obj, VALUE, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
REAL(DFP), INTENT(IN) :: VALUE(:)
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set3
END INTERFACE
- Set all values of scalarfield using a vector of reals.
- The size of value should be same as the size of scalar field.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set3, setting all the values using a vector
PROGRAM main
USE AbstractField_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE GlobalData
USE ScalarField_Class
USE FEDOF_Class
USE Display_Method
USE ApproxUtility
USE Test_Method
USE ReallocateUtility
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
IMPLICIT NONE
CHARACTER(*), PARAMETER :: tomlFileName = "./Set1.toml", &
myName = "main", &
modName = "_ImportFromToml_test_1.F90"
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(FEDOF_) :: fedof, geofedof
TYPE(ScalarField_) :: u
CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL dom%ImportFromToml(fileName=tomlFileName, tomlName="domain")
CALL u%ImportFromToml(fedof=fedof, geofedof=geofedof, dom=dom, &
fileName=tomlFileName, tomlName="u")
mesh => dom%GetMeshPointer()
CALL test1
CALL test2
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
CHARACTER(*), PARAMETER :: testName = "test1()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans
LOGICAL(LGT) :: isok
REAL(DFP), ALLOCATABLE :: realVec(:)
tNodes = u%SIZE()
CALL Reallocate(realVec, tNodes)
CALL RANDOM_NUMBER(realVec)
CALL u%Set(VALUE=realVec)
DO ii = 1, tNodes
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
areal = realVec(ii)
isok = areal.APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test2
!----------------------------------------------------------------------------
SUBROUTINE test2
CHARACTER(*), PARAMETER :: testName = "test2()"
INTEGER(I4B) :: ii, tNodes
REAL(DFP) :: areal, ans, scale
LOGICAL(LGT) :: isok
REAL(DFP), ALLOCATABLE :: realVec(:)
tNodes = u%SIZE()
CALL Reallocate(realVec, tNodes)
CALL RANDOM_NUMBER(realVec)
CALL u%Set(VALUE=0.0_DFP)
scale = 2.0_DFP
CALL u%Set(VALUE=realVec, scale=scale, addContribution=.TRUE.)
DO ii = 1, tNodes
CALL u%Get(globalNode=ii, islocal=.TRUE., VALUE=ans)
areal = realVec(ii) * scale
isok = areal.APPROXEQ.ans
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test2
END PROGRAM main
Interface 4
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set4(obj, globalNode, VALUE, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: globalNode(:)
REAL(DFP), INTENT(IN) :: VALUE
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set4
END INTERFACE
Set multiple values to a scalar value.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set4 and Set5, setting all the values using index
PROGRAM main
USE AbstractField_Class
USE FEDomain_Class
USE AbstractMesh_Class
USE GlobalData
USE ScalarField_Class
USE FEDOF_Class
USE Display_Method
USE ApproxUtility
USE Test_Method
USE ReallocateUtility
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
IMPLICIT NONE
CHARACTER(*), PARAMETER :: tomlFileName = "./Set1.toml", &
myName = "main", &
modName = "_ImportFromToml_test_1.F90"
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(FEDOF_) :: fedof, geofedof
TYPE(ScalarField_) :: u
CALL e%SetQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL dom%ImportFromToml(fileName=tomlFileName, tomlName="domain")
CALL u%ImportFromToml(fedof=fedof, geofedof=geofedof, dom=dom, &
fileName=tomlFileName, tomlName="u")
mesh => dom%GetMeshPointer()
CALL test1
CALL test2
CONTAINS
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test1
CHARACTER(*), PARAMETER :: testName = "test1()"
INTEGER(I4B) :: iel, ii, tNodes, tElements, tcon, maxCon
LOGICAL(LGT) :: isok
REAL(DFP), ALLOCATABLE :: realVec(:), ans(:)
INTEGER(I4B), ALLOCATABLE :: con(:)
tNodes = u%SIZE()
tElements = mesh%GetTotalElements()
maxCon = fedof%GetMaxTotalConnectivity()
CALL Reallocate(con, maxCon)
CALL Reallocate(ans, maxCon)
CALL Reallocate(realVec, maxCon)
CALL RANDOM_NUMBER(realVec)
DO iel = 1, tElements
CALL fedof%GetConnectivity_(ans=con, tsize=tcon, opt="A", &
globalElement=iel, islocal=.TRUE.)
CALL RANDOM_NUMBER(realVec(1:tcon))
CALL u%Set(VALUE=realVec(1:tcon), globalNode=con(1:tcon), &
islocal=.TRUE.)
CALL u%Get(globalNode=con(1:tcon), islocal=.TRUE., &
VALUE=ans, tsize=tcon)
DO ii = 1, tcon
isok = realVec(ii) .APPROXEQ.ans(ii)
IF (.NOT. isok) EXIT
END DO
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test1
!----------------------------------------------------------------------------
! test1
!----------------------------------------------------------------------------
SUBROUTINE test2
CHARACTER(*), PARAMETER :: testName = "test2()"
INTEGER(I4B) :: iel, ii, tNodes, tElements, tcon, maxCon
LOGICAL(LGT) :: isok
REAL(DFP) :: areal
REAL(DFP), ALLOCATABLE :: realVec(:), ans(:)
INTEGER(I4B), ALLOCATABLE :: con(:)
tNodes = u%SIZE()
tElements = mesh%GetTotalElements()
maxCon = fedof%GetMaxTotalConnectivity()
CALL Reallocate(con, maxCon)
CALL Reallocate(ans, maxCon)
DO iel = 1, tElements
CALL fedof%GetConnectivity_(ans=con, tsize=tcon, opt="A", &
globalElement=iel, islocal=.TRUE.)
CALL RANDOM_NUMBER(areal)
CALL u%Set(VALUE=areal, globalNode=con(1:tcon), &
islocal=.TRUE.)
CALL u%Get(globalNode=con(1:tcon), islocal=.TRUE., &
VALUE=ans, tsize=tcon)
DO ii = 1, tcon
isok = areal.APPROXEQ.ans(ii)
IF (.NOT. isok) EXIT
END DO
IF (.NOT. isok) THEN
CALL Display("Test Failed in "//testName//": at node "//ToString(ii))
RETURN
END IF
END DO
CALL OK(.TRUE., testName)
END SUBROUTINE test2
!----------------------------------------------------------------------------
!
!----------------------------------------------------------------------------
END PROGRAM main
Interface 5
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set5(obj, globalNode, VALUE, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: globalNode(:)
REAL(DFP), INTENT(IN) :: VALUE(:)
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set5
END INTERFACE
Set multiple values of scalar field.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Testing Set5
PROGRAM main
USE FEDomain_Class
USE HDF5File_Class
USE AbstractMesh_Class
USE AbstractField_Class, ONLY: TypeField
USE ScalarField_Class
USE FPL, ONLY: FPL_Init, FPL_FINALIZE, ParameterList_
USE GlobalData
USE Test_Method
USE FEDOF_Class
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ApproxUtility
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(HDF5File_) :: meshfile
TYPE(ParameterList_) :: param
CHARACTER(LEN=*), PARAMETER :: engine = "NATIVE_SERIAL"
CHARACTER(*), PARAMETER :: meshfilename = &
"../../Mesh/examples/meshdata/small_mesh.h5"
INTEGER(I4B), PARAMETER :: nsd = 2
CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL FPL_Init()
CALL param%Initiate()
!> start creating domain
CALL meshfile%Initiate(filename=meshfilename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(hdf5=meshfile, group="")
!> end creating domain
mesh => dom%GetMeshPointer(dim=nsd)
BLOCK
INTEGER(I4B), PARAMETER :: order = 1
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Lagrange"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj
REAL(DFP) :: found(100), want(100), VALUE(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set5 "
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE., &
addContribution=.TRUE., scale=2.0_DFP)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = 3 * VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE., &
addContribution=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = 4 * VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
BLOCK
INTEGER(I4B), PARAMETER :: order = 2
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Heirarchical"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj
REAL(DFP) :: found(100), want(100), VALUE(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set5 "
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
mesh => NULL()
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CALL param%DEALLOCATE()
CALL FPL_FINALIZE()
END PROGRAM main
Interface 6
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set6(obj, istart, iend, stride, VALUE, &
& scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: istart
INTEGER(I4B), INTENT(IN) :: iend
INTEGER(I4B), INTENT(IN) :: stride
REAL(DFP), INTENT(IN) :: VALUE
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set6
END INTERFACE
Set multiple values by using triplets istart:iend:stride
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set5
PROGRAM main
USE FEDomain_Class
USE HDF5File_Class
USE AbstractMesh_Class
USE AbstractField_Class, ONLY: TypeField
USE ScalarField_Class
USE FPL, ONLY: FPL_Init, FPL_FINALIZE, ParameterList_
USE GlobalData
USE Test_Method
USE FEDOF_Class
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ApproxUtility
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(HDF5File_) :: meshfile
TYPE(ParameterList_) :: param
CHARACTER(LEN=*), PARAMETER :: engine = "NATIVE_SERIAL"
CHARACTER(*), PARAMETER :: meshfilename = &
"../../Mesh/examples/meshdata/small_mesh.h5"
INTEGER(I4B), PARAMETER :: nsd = 2
CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL FPL_Init()
CALL param%Initiate()
!> start creating domain
CALL meshfile%Initiate(filename=meshfilename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(hdf5=meshfile, group="")
!> end creating domain
mesh => dom%GetMeshPointer(dim=nsd)
BLOCK
INTEGER(I4B), PARAMETER :: order = 1
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Lagrange"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj
REAL(DFP) :: found(100), want(100), VALUE(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set5 "
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE., &
addContribution=.TRUE., scale=2.0_DFP)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = 3 * VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE., &
addContribution=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = 4 * VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
BLOCK
INTEGER(I4B), PARAMETER :: order = 2
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Heirarchical"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj
REAL(DFP) :: found(100), want(100), VALUE(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set5 "
VALUE(1:3) = [10, 20, 30]
CALL obj%Set(VALUE=VALUE(1:3), globalNode=localNode, islocal=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = VALUE(1:3)
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
mesh => NULL()
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CALL param%DEALLOCATE()
CALL FPL_FINALIZE()
END PROGRAM main
Interface 7
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set7(obj, istart, iend, stride, VALUE, &
& scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: istart
INTEGER(I4B), INTENT(IN) :: iend
INTEGER(I4B), INTENT(IN) :: stride
REAL(DFP), INTENT(IN) :: VALUE(:)
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set7
END INTERFACE
Set multiple values using triplets.
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set7
PROGRAM main
USE FEDomain_Class
USE HDF5File_Class
USE AbstractMesh_Class
USE AbstractField_Class, ONLY: TypeField
USE ScalarField_Class
USE FPL, ONLY: FPL_Init, FPL_FINALIZE, ParameterList_
USE GlobalData
USE Test_Method
USE FEDOF_Class
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ApproxUtility
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(HDF5File_) :: meshfile
TYPE(ParameterList_) :: param
CHARACTER(LEN=*), PARAMETER :: engine = "NATIVE_SERIAL"
CHARACTER(*), PARAMETER :: meshfilename = &
"../../Mesh/examples/meshdata/small_mesh.h5"
INTEGER(I4B), PARAMETER :: nsd = 2
CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL FPL_Init()
CALL param%Initiate()
!> start creating domain
CALL meshfile%Initiate(filename=meshfilename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(hdf5=meshfile, group="")
!> end creating domain
mesh => dom%GetMeshPointer(dim=nsd)
BLOCK
INTEGER(I4B), PARAMETER :: order = 1
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Lagrange"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj, VALUE
REAL(DFP) :: found(100), want(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
CALL VALUE%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set7 "
DO ii = 1, SIZE(localNode)
CALL VALUE%Set(VALUE=REAL(localNode(ii), kind=DFP), &
globalNode=localNode(ii), islocal=.TRUE.)
END DO
CALL obj%Set(VALUE=VALUE)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = localNode
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%Set(0.0_DFP)
obj = VALUE
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = localNode
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
mesh => NULL()
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CALL param%DEALLOCATE()
CALL FPL_FINALIZE()
END PROGRAM main
Interface 8
- ܀ Interface
- ️܀ See example
- ↢
INTERFACE
MODULE SUBROUTINE set8(obj, obj2)
CLASS(ScalarField_), INTENT(INOUT) :: obj
CLASS(ScalarField_), INTENT(IN) :: obj2
END SUBROUTINE set8
END INTERFACE
!> author: Vikas Sharma, Ph. D.
! date: 2024-06-05
! summary: Set8
PROGRAM main
USE FEDomain_Class
USE HDF5File_Class
USE AbstractMesh_Class
USE AbstractField_Class, ONLY: TypeField
USE ScalarField_Class
USE FPL, ONLY: FPL_Init, FPL_FINALIZE, ParameterList_
USE GlobalData
USE Test_Method
USE FEDOF_Class
USE ExceptionHandler_Class, ONLY: e, EXCEPTION_INFORMATION
USE ApproxUtility
TYPE(FEDomain_) :: dom
CLASS(AbstractMesh_), POINTER :: mesh
TYPE(HDF5File_) :: meshfile
TYPE(ParameterList_) :: param
CHARACTER(LEN=*), PARAMETER :: engine = "NATIVE_SERIAL"
CHARACTER(*), PARAMETER :: meshfilename = &
"../../Mesh/examples/meshdata/small_mesh.h5"
INTEGER(I4B), PARAMETER :: nsd = 2
CALL e%setQuietMode(EXCEPTION_INFORMATION, .TRUE.)
CALL FPL_Init()
CALL param%Initiate()
!> start creating domain
CALL meshfile%Initiate(filename=meshfilename, mode="READ")
CALL meshfile%OPEN()
CALL dom%Initiate(hdf5=meshfile, group="")
!> end creating domain
mesh => dom%GetMeshPointer(dim=nsd)
BLOCK
INTEGER(I4B), PARAMETER :: order = 1
CHARACTER(*), PARAMETER :: baseContinuity = "H1"
CHARACTER(*), PARAMETER :: baseInterpolation = "Lagrange"
TYPE(FEDOF_) :: fedof
TYPE(ScalarField_) :: obj, VALUE
REAL(DFP) :: found(100), want(100), tol
INTEGER(I4B) :: tsize, localNode(3)
CHARACTER(:), ALLOCATABLE :: msg
LOGICAL(LGT) :: isok
CALL fedof%Initiate(baseContinuity=baseContinuity, &
baseInterpolation=baseInterpolation, order=order, mesh=mesh)
CALL SetScalarFieldParam(param=param, &
fieldType=TypeField%normal, &
name="U", &
engine=engine)
CALL obj%Initiate(param, fedof)
CALL VALUE%Initiate(param, fedof)
localNode = [1, 3, 5]
msg = "Set8 "
DO ii = 1, SIZE(localNode)
CALL VALUE%Set(VALUE=REAL(localNode(ii), kind=DFP), &
globalNode=localNode(ii), islocal=.TRUE.)
END DO
CALL obj%Set(VALUE=VALUE)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = localNode
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%Set(VALUE=VALUE, scale=2.0_DFP, addContribution=.TRUE.)
CALL obj%Get(VALUE=found, tsize=tsize)
want = 0.0_DFP
want(localNode) = 3 * localNode
tol = 1.0E-5
isok = ALL(SOFTEQ(found(1:tsize), want(1:tsize), tol))
CALL OK(isok, msg)
CALL obj%DEALLOCATE()
END BLOCK
mesh => NULL()
CALL dom%DEALLOCATE()
CALL meshfile%DEALLOCATE()
CALL param%DEALLOCATE()
CALL FPL_FINALIZE()
END PROGRAM main
Interface 9
INTERFACE
MODULE SUBROUTINE set9(obj, globalNode, VALUE, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: globalNode(:)
TYPE(FEVariable_), INTENT(IN) :: VALUE
!! Scalar, Nodal, FEVariable (Space or Constant)
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE set9
END INTERFACE
Interface 10
INTERFACE
MODULE SUBROUTINE set10(obj, obj2, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
CLASS(ScalarField_), INTENT(IN) :: obj2
REAL(DFP), INTENT(IN) :: scale
LOGICAL(LGT), INTENT(IN) :: addContribution
END SUBROUTINE set10
END INTERFACE
Interface 11
INTERFACE
MODULE SUBROUTINE Set11(obj, ivar, idof, VALUE, ivar_value, &
& idof_value, scale, addContribution)
CLASS(ScalarField_), INTENT(INOUT) :: obj
INTEGER(I4B), INTENT(IN) :: ivar
INTEGER(I4B), INTENT(IN) :: idof
CLASS(AbstractNodeField_), INTENT(IN) :: VALUE
INTEGER(I4B), INTENT(IN) :: ivar_value
INTEGER(I4B), INTENT(IN) :: idof_value
REAL(DFP), OPTIONAL, INTENT(IN) :: scale
LOGICAL(LGT), OPTIONAL, INTENT(IN) :: addContribution
END SUBROUTINE Set11
END INTERFACE