DAMASK_EICMD/src/mesh/FEM_utilities.f90

205 lines
8.5 KiB
Fortran
Raw Normal View History

2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Utilities used by the FEM solver
!--------------------------------------------------------------------------------------------------
module FEM_utilities
2018-08-30 16:08:44 +05:30
#include <petsc/finclude/petscdmplex.h>
2018-09-21 11:49:36 +05:30
#include <petsc/finclude/petscdmda.h>
#include <petsc/finclude/petscis.h>
use PETScDMplex
use PETScDMDA
use PETScIS
#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY)
use MPI_f08
#endif
2019-06-11 13:18:07 +05:30
use prec
2020-08-15 19:32:10 +05:30
use config
2019-06-11 13:18:07 +05:30
use math
use IO
2020-03-20 19:38:07 +05:30
use discretization_mesh
use homogenization
use FEM_quadrature
2019-06-11 13:18:07 +05:30
implicit none
private
2018-08-17 03:44:25 +05:30
2021-07-22 11:18:01 +05:30
logical, public :: cutBack = .false. !< cut back of BVP solver in case convergence is not achieved or a material point is terminally ill
real(pReal), public, protected :: wgt !< weighting factor 1/Nelems
2019-05-12 13:35:40 +05:30
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! field labels information
character(len=*), parameter, public :: &
FIELD_MECH_label = 'mechanical'
2020-03-17 12:47:14 +05:30
enum, bind(c); enumerator :: &
FIELD_UNDEFINED_ID, &
FIELD_MECH_ID
2019-06-11 13:18:07 +05:30
end enum
2020-03-17 12:47:14 +05:30
enum, bind(c); enumerator :: &
COMPONENT_UNDEFINED_ID, &
COMPONENT_MECH_X_ID, &
COMPONENT_MECH_Y_ID, &
COMPONENT_MECH_Z_ID
2019-06-11 13:18:07 +05:30
end enum
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! derived types
2019-06-11 13:18:07 +05:30
type, public :: tSolutionState !< return type of solution from FEM solver variants
logical :: converged = .true.
logical :: stagConverged = .true.
PetscInt :: iterationsNeeded = 0_pPETSCINT
2019-06-11 13:18:07 +05:30
end type tSolutionState
2019-06-11 13:18:07 +05:30
type, public :: tComponentBC
integer(kind(COMPONENT_UNDEFINED_ID)) :: ID
real(pReal), allocatable, dimension(:) :: Value
logical, allocatable, dimension(:) :: Mask
2019-06-11 13:18:07 +05:30
end type tComponentBC
2019-06-11 13:18:07 +05:30
type, public :: tFieldBC
integer(kind(FIELD_UNDEFINED_ID)) :: ID
integer :: nComponents = 0
type(tComponentBC), allocatable, dimension(:) :: componentBC
2019-06-11 13:18:07 +05:30
end type tFieldBC
2019-06-11 13:18:07 +05:30
public :: &
2020-06-17 21:32:22 +05:30
FEM_utilities_init, &
2019-06-11 13:18:07 +05:30
utilities_constitutiveResponse, &
utilities_projectBCValues, &
FIELD_MECH_ID, &
2020-01-21 12:16:32 +05:30
COMPONENT_UNDEFINED_ID, &
2019-06-11 13:18:07 +05:30
COMPONENT_MECH_X_ID, &
COMPONENT_MECH_Y_ID, &
COMPONENT_MECH_Z_ID
2018-08-17 03:44:25 +05:30
contains
2018-08-17 03:44:25 +05:30
2020-06-29 18:39:13 +05:30
!ToDo: use functions in variable call
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields, sets debug flags
!--------------------------------------------------------------------------------------------------
2020-06-17 21:32:22 +05:30
subroutine FEM_utilities_init
2020-01-26 16:54:35 +05:30
character(len=pStringLen) :: petsc_optionsOrder
2020-06-16 22:45:01 +05:30
class(tNode), pointer :: &
2020-06-18 02:30:03 +05:30
num_mesh, &
2020-06-26 23:42:05 +05:30
debug_mesh ! pointer to mesh debug options
integer :: &
p_s, & !< order of shape functions
p_i !< integration order (quadrature rule)
2020-07-03 20:15:11 +05:30
character(len=*), parameter :: &
2020-06-18 21:44:53 +05:30
PETSCDEBUG = ' -snes_view -snes_monitor '
PetscErrorCode :: err_PETSc
2021-07-22 11:18:01 +05:30
logical :: debugPETSc !< use some in debug defined options for more verbose PETSc solution
2018-08-17 03:44:25 +05:30
print'(/,1x,a)', '<<<+- FEM_utilities init -+>>>'
2020-09-13 14:47:49 +05:30
num_mesh => config_numerics%get('mesh',defaultVal=emptyDict)
p_s = num_mesh%get_asInt('p_s',defaultVal = 2)
p_i = num_mesh%get_asInt('p_i',defaultVal = p_s)
2022-01-13 12:13:53 +05:30
if (p_s < 1 .or. p_s > size(FEM_nQuadrature,2)) &
call IO_error(821,ext_msg='shape function order (p_s) out of bounds')
2022-01-13 12:13:53 +05:30
if (p_i < max(1,p_s-1) .or. p_i > p_s) &
call IO_error(821,ext_msg='integration order (p_i) out of bounds')
2020-06-16 22:45:01 +05:30
2020-09-13 14:47:49 +05:30
debug_mesh => config_debug%get('mesh',defaultVal=emptyList)
2021-03-26 17:04:16 +05:30
debugPETSc = debug_mesh%contains('PETSc')
2020-07-03 20:15:11 +05:30
if(debugPETSc) print'(3(/,1x,a),/)', &
'Initializing PETSc with debug options: ', &
2019-06-11 13:18:07 +05:30
trim(PETScDebug), &
'add more using the "PETSc_options" keyword in numerics.yaml'
2020-09-22 16:39:12 +05:30
flush(IO_STDOUT)
call PetscOptionsClear(PETSC_NULL_OPTIONS,err_PETSc)
CHKERRQ(err_PETSc)
if(debugPETSc) call PetscOptionsInsertString(PETSC_NULL_OPTIONS,trim(PETSCDEBUG),err_PETSc)
CHKERRQ(err_PETSc)
2021-02-09 03:51:53 +05:30
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-mechanical_snes_type newtonls &
&-mechanical_snes_linesearch_type cp -mechanical_snes_ksp_ew &
&-mechanical_snes_ksp_ew_rtol0 0.01 -mechanical_snes_ksp_ew_rtolmax 0.01 &
&-mechanical_ksp_type fgmres -mechanical_ksp_max_it 25', err_PETSc)
CHKERRQ(err_PETSc)
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_mesh%get_asString('PETSc_options',defaultVal=''),err_PETSc)
CHKERRQ(err_PETSc)
write(petsc_optionsOrder,'(a,i0)') '-mechFE_petscspace_degree ', p_s
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,trim(petsc_optionsOrder),err_PETSc)
CHKERRQ(err_PETSc)
2019-06-11 13:18:07 +05:30
wgt = 1.0/real(mesh_maxNips*mesh_NcpElemsGlobal,pReal)
2018-08-17 03:44:25 +05:30
2020-06-17 21:32:22 +05:30
end subroutine FEM_utilities_init
2018-08-17 03:44:25 +05:30
2020-03-17 04:40:23 +05:30
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates constitutive response
!--------------------------------------------------------------------------------------------------
subroutine utilities_constitutiveResponse(timeinc,P_av,forwardData)
2019-06-11 13:18:07 +05:30
real(pReal), intent(in) :: timeinc !< loading time
logical, intent(in) :: forwardData !< age results
real(pReal),intent(out), dimension(3,3) :: P_av !< average PK stress
integer(MPI_INTEGER_KIND) :: err_MPI
2018-08-17 03:44:25 +05:30
print'(/,1x,a)', '... evaluating constitutive response ......................................'
2018-08-17 03:44:25 +05:30
call homogenization_mechanical_response(timeinc,1,mesh_maxNips*mesh_NcpElems) ! calculate P field
2021-07-17 00:13:08 +05:30
if (.not. terminallyIll) &
2021-07-17 15:20:21 +05:30
call homogenization_mechanical_response2(timeinc,[1,mesh_maxNips],[1,mesh_NcpElems])
2021-07-17 00:13:08 +05:30
cutBack = .false.
P_av = sum(homogenization_P,dim=3) * wgt
call MPI_Allreduce(MPI_IN_PLACE,P_av,9_MPI_INTEGER_KIND,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) error stop 'MPI error'
2018-08-17 03:44:25 +05:30
end subroutine utilities_constitutiveResponse
2019-03-09 03:53:07 +05:30
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Project BC values to local vector
!--------------------------------------------------------------------------------------------------
subroutine utilities_projectBCValues(localVec,section,field,comp,bcPointsIS,BCValue,BCDotValue,timeinc)
2019-06-11 13:18:07 +05:30
Vec :: localVec
PetscInt :: field, comp, nBcPoints, point, dof, numDof, numComp, offset
PetscSection :: section
IS :: bcPointsIS
PetscInt, pointer :: bcPoints(:)
PetscScalar, pointer :: localArray(:)
PetscScalar :: BCValue,BCDotValue,timeinc
PetscErrorCode :: err_PETSc
2019-06-11 13:18:07 +05:30
2021-11-12 04:29:48 +05:30
call PetscSectionGetFieldComponents(section,field,numComp,err_PETSc)
CHKERRQ(err_PETSc)
call ISGetSize(bcPointsIS,nBcPoints,err_PETSc)
CHKERRQ(err_PETSc)
if (nBcPoints > 0) call ISGetIndicesF90(bcPointsIS,bcPoints,err_PETSc)
call VecGetArrayF90(localVec,localArray,err_PETSc); CHKERRQ(err_PETSc)
2019-06-11 13:18:07 +05:30
do point = 1, nBcPoints
call PetscSectionGetFieldDof(section,bcPoints(point),field,numDof,err_PETSc)
CHKERRQ(err_PETSc)
call PetscSectionGetFieldOffset(section,bcPoints(point),field,offset,err_PETSc)
CHKERRQ(err_PETSc)
2019-06-11 13:18:07 +05:30
do dof = offset+comp+1, offset+numDof, numComp
localArray(dof) = localArray(dof) + BCValue + BCDotValue*timeinc
2021-11-12 04:29:48 +05:30
end do
end do
call VecRestoreArrayF90(localVec,localArray,err_PETSc); CHKERRQ(err_PETSc)
call VecAssemblyBegin(localVec, err_PETSc); CHKERRQ(err_PETSc)
call VecAssemblyEnd (localVec, err_PETSc); CHKERRQ(err_PETSc)
if (nBcPoints > 0) call ISRestoreIndicesF90(bcPointsIS,bcPoints,err_PETSc)
CHKERRQ(err_PETSc)
2018-08-17 03:44:25 +05:30
end subroutine utilities_projectBCValues
end module FEM_utilities