LinearStaticCDR_test_6
This example demonstrates how to USE addMaterial
to add NEW material to LinearStaticCDR_
kernel.
PROGRAM main
USE easifemBase
USE easifemClasses
USE easifemMaterials
USE easifemKernels
USE LinearStaticCDR_Class
IMPLICIT NONE
TYPE(LinearStaticCDR_) :: obj
TYPE(ParameterList_) :: param
TYPE(Domain_) :: dom
CHARACTER( LEN = * ), PARAMETER :: domainFileName="./mesh.h5"
TYPE( HDF5File_ ) :: domainFile
TYPE( MeshSelection_ ) :: region
INTEGER(I4B), PARAMETER :: tMaterials = 2, tDirichletBC=4, &
& solverName = LIS_GMRES, preconditionOption=LEFT_PRECONDITION, &
& convergenceIn = convergenceInRes, &
& convergenceType=relativeConvergence, &
& maxIter = 100, &
& KrylovSubspaceSize=20
!! Initiate an instance of [[ParameterList_]]
CALL FPL_INIT(); CALL param%Initiate()
!! Set the PARAMETER for [[LinearStaticCDR_]]
CALL SetLinearStaticCDRParam( param=param, &
& engine="NATIVE_SERIAL", &
& isConservative=.TRUE., &
& coordinateSystem=KERNEL_2D, &
& tMaterials=tMaterials, &
& tDirichletBC=tDirichletBC, &
& domainFile=domainFileName, &
& baseContinuity="H1", &
& baseInterpolation="LagrangeInterpolation", &
& quadratureType="GaussLegendre" )
!! Set the PARAMETER for [[LinSolver_]].
CALL SetLinSolverParam( &
& param=param, &
& solverName=solverName,&
& preconditionOption=preconditionOption, &
& convergenceIn=convergenceIn, &
& convergenceType=convergenceType, &
& maxIter=maxIter, &
& relativeToRHS=.TRUE., &
& KrylovSubspaceSize=KrylovSubspaceSize, &
& rtol=1.0D-10, &
& atol=1.0D-10 )
!! Initiates computation domain.
CALL domainFile%Initiate(filename=domainFileName, mode="READ")
CALL domainFile%Open()
CALL dom%Initiate( domainFile, '' )
CALL domainFile%Deallocate()
!! Initiate an instace of [[LinearStaticCDR_]] kernel
CALL obj%Initiate( param=param, dom=dom )
!! Add another material and domain region in the [[LinearStaticCDR_]] kernel.
CALL region%Initiate( isSelectionByMeshID=.TRUE. )
CALL region%Add( dim=2, meshID=[1] )
CALL SetSolidMaterialParam( param=param, name="SolidMaterial", &
& massDensity=7000.0_DFP, diffusivity=11.72D-6 )
CALL obj%AddMaterial( materialNo=1, materialName="SolidMaterial", &
& param=param, region=region )
CALL region%Deallocate()
!! Add another material and domain region in the [[LinearStaticCDR_]] kernel.
CALL region%Initiate( isSelectionByMeshID=.TRUE. )
CALL region%Add( dim=2, meshID=[2] )
CALL SetSolidMaterialParam( param=param, name="SolidMaterial", &
& massDensity=7000.0_DFP, diffusivity=11.72D-6 )
CALL obj%AddMaterial( materialNo=2, materialName="SolidMaterial", &
& param=param, region=region )
CALL region%Deallocate()
!! Display the content of kernel.
CALL obj%Display( "LinearStaticCDR :: ")
CALL obj%Deallocate( )
CALL dom%Deallocate()
CALL param%Deallocate()
CALL FPL_FINALIZE()
END PROGRAM main