2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
2019-03-12 10:36:59 +05:30
|
|
|
|
!> @brief Grid solver for mechanics: Spectral basic
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-12 10:36:59 +05:30
|
|
|
|
module grid_mech_spectral_basic
|
2018-05-17 15:34:21 +05:30
|
|
|
|
#include <petsc/finclude/petscsnes.h>
|
|
|
|
|
#include <petsc/finclude/petscdmda.h>
|
2019-03-12 16:06:18 +05:30
|
|
|
|
use PETScdmda
|
|
|
|
|
use PETScsnes
|
2019-06-10 00:50:38 +05:30
|
|
|
|
|
|
|
|
|
use prec
|
|
|
|
|
use DAMASK_interface
|
|
|
|
|
use HDF5_utilities
|
|
|
|
|
use math
|
|
|
|
|
use spectral_utilities
|
|
|
|
|
use FEsolving
|
|
|
|
|
use config
|
|
|
|
|
use numerics
|
|
|
|
|
use homogenization
|
2020-03-20 11:28:55 +05:30
|
|
|
|
use discretization_grid
|
2019-06-10 00:50:38 +05:30
|
|
|
|
use debug
|
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
implicit none
|
|
|
|
|
private
|
2019-03-23 15:16:56 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! derived types
|
2019-03-12 16:06:18 +05:30
|
|
|
|
type(tSolutionParams), private :: params
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-26 12:06:55 +05:30
|
|
|
|
type, private :: tNumerics
|
2019-04-15 19:23:46 +05:30
|
|
|
|
logical :: update_gamma !< update gamma operator with current stiffness
|
2019-03-26 12:06:55 +05:30
|
|
|
|
end type tNumerics
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-26 12:06:55 +05:30
|
|
|
|
type(tNumerics) :: num ! numerics parameters. Better name?
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! PETSc data
|
2019-03-12 16:06:18 +05:30
|
|
|
|
DM, private :: da
|
|
|
|
|
SNES, private :: snes
|
|
|
|
|
Vec, private :: solution_vec
|
2012-11-12 19:44:39 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! common pointwise data
|
2019-03-22 20:32:00 +05:30
|
|
|
|
real(pReal), private, dimension(:,:,:,:,:), allocatable :: &
|
2019-04-15 19:23:46 +05:30
|
|
|
|
F_lastInc, & !< field of previous compatible deformation gradients
|
|
|
|
|
Fdot !< field of assumed rate of compatible deformation gradient
|
2012-11-12 19:44:39 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! stress, stiffness and compliance average etc.
|
2019-03-12 16:06:18 +05:30
|
|
|
|
real(pReal), private, dimension(3,3) :: &
|
|
|
|
|
F_aimDot = 0.0_pReal, & !< assumed rate of average deformation gradient
|
|
|
|
|
F_aim = math_I3, & !< current prescribed deformation gradient
|
|
|
|
|
F_aim_lastInc = math_I3, & !< previous average deformation gradient
|
|
|
|
|
P_av = 0.0_pReal !< average 1st Piola--Kirchhoff stress
|
2019-03-22 20:32:00 +05:30
|
|
|
|
|
2020-02-25 23:23:12 +05:30
|
|
|
|
character(len=pStringLen), private :: incInfo !< time and increment information
|
2019-03-12 16:06:18 +05:30
|
|
|
|
real(pReal), private, dimension(3,3,3,3) :: &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
C_volAvg = 0.0_pReal, & !< current volume average stiffness
|
2019-03-12 16:06:18 +05:30
|
|
|
|
C_volAvgLastInc = 0.0_pReal, & !< previous volume average stiffness
|
|
|
|
|
C_minMaxAvg = 0.0_pReal, & !< current (min+max)/2 stiffness
|
|
|
|
|
C_minMaxAvgLastInc = 0.0_pReal, & !< previous (min+max)/2 stiffness
|
|
|
|
|
S = 0.0_pReal !< current compliance (filled up with zeros)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
real(pReal), private :: &
|
|
|
|
|
err_BC, & !< deviation from stress BC
|
|
|
|
|
err_div !< RMS of div of P
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
integer, private :: &
|
|
|
|
|
totalIter = 0 !< total iteration in current increment
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
public :: &
|
|
|
|
|
grid_mech_spectral_basic_init, &
|
|
|
|
|
grid_mech_spectral_basic_solution, &
|
2019-10-24 17:16:36 +05:30
|
|
|
|
grid_mech_spectral_basic_forward, &
|
2019-10-28 18:06:36 +05:30
|
|
|
|
grid_mech_spectral_basic_updateCoords, &
|
2019-10-24 17:16:36 +05:30
|
|
|
|
grid_mech_spectral_basic_restartWrite
|
|
|
|
|
|
2013-01-10 03:49:32 +05:30
|
|
|
|
contains
|
2012-11-09 03:03:58 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-02-16 20:06:18 +05:30
|
|
|
|
!> @brief allocates all necessary fields and fills them with data, potentially from restart info
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-12 10:36:59 +05:30
|
|
|
|
subroutine grid_mech_spectral_basic_init
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
real(pReal), dimension(3,3,grid(1),grid(2),grid3) :: P
|
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
|
temp33_Real = 0.0_pReal
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
PetscErrorCode :: ierr
|
2019-03-22 20:32:00 +05:30
|
|
|
|
PetscScalar, pointer, dimension(:,:,:,:) :: &
|
2019-04-10 16:53:57 +05:30
|
|
|
|
F ! pointer to solution data
|
2020-02-25 23:23:12 +05:30
|
|
|
|
PetscInt, dimension(worldsize) :: localK
|
2019-12-07 19:54:45 +05:30
|
|
|
|
integer(HID_T) :: fileHandle, groupHandle
|
|
|
|
|
integer :: fileUnit
|
|
|
|
|
character(len=pStringLen) :: fileName
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(6,'(/,a)') ' <<<+- grid_mech_spectral_basic init -+>>>'; flush(6)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
write(6,'(/,a)') ' Eisenlohr et al., International Journal of Plasticity 46:37–53, 2013'
|
|
|
|
|
write(6,'(a)') ' https://doi.org/10.1016/j.ijplas.2012.09.012'
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
write(6,'(/,a)') ' Shanthraj et al., International Journal of Plasticity 66:31–45, 2015'
|
|
|
|
|
write(6,'(a)') ' https://doi.org/10.1016/j.ijplas.2014.02.006'
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-26 12:06:55 +05:30
|
|
|
|
num%update_gamma = config_numerics%getInt('update_gamma',defaultVal=0) > 0
|
2013-12-18 14:39:32 +05:30
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! set default and user defined options for PETSc
|
2019-03-12 16:06:18 +05:30
|
|
|
|
call PETScOptionsInsertString(PETSC_NULL_OPTIONS,'-mech_snes_type ngmres',ierr)
|
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
|
call PETScOptionsInsertString(PETSC_NULL_OPTIONS,trim(petsc_options),ierr)
|
|
|
|
|
CHKERRQ(ierr)
|
2019-03-12 10:36:59 +05:30
|
|
|
|
|
2012-11-12 19:44:39 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! allocate global fields
|
2019-03-12 16:06:18 +05:30
|
|
|
|
allocate (F_lastInc(3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
|
|
|
|
|
allocate (Fdot (3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2012-11-12 19:44:39 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! initialize solver specific parts of PETSc
|
2019-03-12 16:06:18 +05:30
|
|
|
|
call SNESCreate(PETSC_COMM_WORLD,snes,ierr); CHKERRQ(ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
call SNESSetOptionsPrefix(snes,'mech_',ierr);CHKERRQ(ierr)
|
2019-03-12 16:06:18 +05:30
|
|
|
|
localK = 0
|
|
|
|
|
localK(worldrank+1) = grid3
|
|
|
|
|
call MPI_Allreduce(MPI_IN_PLACE,localK,worldsize,MPI_INTEGER,MPI_SUM,PETSC_COMM_WORLD,ierr)
|
|
|
|
|
call DMDACreate3d(PETSC_COMM_WORLD, &
|
|
|
|
|
DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, & ! cut off stencil at boundary
|
|
|
|
|
DMDA_STENCIL_BOX, & ! Moore (26) neighborhood around central point
|
|
|
|
|
grid(1),grid(2),grid(3), & ! global grid
|
|
|
|
|
1 , 1, worldsize, &
|
|
|
|
|
9, 0, & ! #dof (F tensor), ghost boundary width (domain overlap)
|
|
|
|
|
[grid(1)],[grid(2)],localK, & ! local grid
|
|
|
|
|
da,ierr) ! handle, error
|
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
|
call SNESSetDM(snes,da,ierr); CHKERRQ(ierr) ! connect snes to da
|
|
|
|
|
call DMsetFromOptions(da,ierr); CHKERRQ(ierr)
|
|
|
|
|
call DMsetUp(da,ierr); CHKERRQ(ierr)
|
|
|
|
|
call DMcreateGlobalVector(da,solution_vec,ierr); CHKERRQ(ierr) ! global solution vector (grid x 9, i.e. every def grad tensor)
|
2019-03-22 20:32:00 +05:30
|
|
|
|
call DMDASNESsetFunctionLocal(da,INSERT_VALUES,formResidual,PETSC_NULL_SNES,ierr) ! residual vector of same shape as solution vector
|
2020-02-25 23:23:12 +05:30
|
|
|
|
CHKERRQ(ierr)
|
2019-04-15 19:46:14 +05:30
|
|
|
|
call SNESsetConvergenceTest(snes,converged,PETSC_NULL_SNES,PETSC_NULL_FUNCTION,ierr) ! specify custom convergence check function "converged"
|
2019-03-12 16:06:18 +05:30
|
|
|
|
CHKERRQ(ierr)
|
|
|
|
|
call SNESsetFromOptions(snes,ierr); CHKERRQ(ierr) ! pull it all together with additional CLI arguments
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
2012-11-12 19:44:39 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-25 23:23:12 +05:30
|
|
|
|
! init fields
|
2019-03-22 20:32:00 +05:30
|
|
|
|
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr) ! places pointer on PETSc data
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
restartRead: if (interface_restartInc > 0) then
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(6,'(/,a,i0,a)') ' reading restart data of increment ', interface_restartInc, ' from file'
|
2019-03-24 14:21:47 +05:30
|
|
|
|
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(fileName,'(a,a,i0,a)') trim(getSolverJobName()),'_',worldrank,'.hdf5'
|
|
|
|
|
fileHandle = HDF5_openFile(fileName)
|
2019-10-28 17:47:05 +05:30
|
|
|
|
groupHandle = HDF5_openGroup(fileHandle,'solver')
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_read(groupHandle,F_aim, 'F_aim')
|
|
|
|
|
call HDF5_read(groupHandle,F_aim_lastInc,'F_aim_lastInc')
|
|
|
|
|
call HDF5_read(groupHandle,F_aimDot, 'F_aimDot')
|
|
|
|
|
call HDF5_read(groupHandle,F, 'F')
|
|
|
|
|
call HDF5_read(groupHandle,F_lastInc, 'F_lastInc')
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-06-15 19:18:47 +05:30
|
|
|
|
elseif (interface_restartInc == 0) then restartRead
|
2019-03-12 16:06:18 +05:30
|
|
|
|
F_lastInc = spread(spread(spread(math_I3,3,grid(1)),4,grid(2)),5,grid3) ! initialize to identity
|
|
|
|
|
F = reshape(F_lastInc,[9,grid(1),grid(2),grid3])
|
2019-06-15 19:18:47 +05:30
|
|
|
|
endif restartRead
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
materialpoint_F0 = reshape(F_lastInc, [3,3,1,product(grid(1:2))*grid3]) ! set starting condition for materialpoint_stressAndItsTangent
|
2019-09-28 03:32:36 +05:30
|
|
|
|
call Utilities_updateCoords(reshape(F,shape(F_lastInc)))
|
2019-03-12 16:06:18 +05:30
|
|
|
|
call Utilities_constitutiveResponse(P,temp33_Real,C_volAvg,C_minMaxAvg, & ! stress field, stress avg, global average of stiffness and (min+max)/2
|
|
|
|
|
reshape(F,shape(F_lastInc)), & ! target F
|
2019-12-03 00:53:50 +05:30
|
|
|
|
0.0_pReal) ! time increment
|
2019-03-22 20:32:00 +05:30
|
|
|
|
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr) ! deassociate pointer
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-06-15 19:18:47 +05:30
|
|
|
|
restartRead2: if (interface_restartInc > 0) then
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(6,'(/,a,i0,a)') ' reading more restart data of increment ', interface_restartInc, ' from file'
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_read(groupHandle,C_volAvg, 'C_volAvg')
|
|
|
|
|
call HDF5_read(groupHandle,C_volAvgLastInc,'C_volAvgLastInc')
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_closeGroup(groupHandle)
|
2019-04-10 16:53:57 +05:30
|
|
|
|
call HDF5_closeFile(fileHandle)
|
|
|
|
|
|
2019-09-23 04:04:05 +05:30
|
|
|
|
call MPI_File_open(PETSC_COMM_WORLD, trim(getSolverJobName())//'.C_ref', &
|
|
|
|
|
MPI_MODE_RDONLY,MPI_INFO_NULL,fileUnit,ierr)
|
|
|
|
|
call MPI_File_read(fileUnit,C_minMaxAvg,81,MPI_DOUBLE,MPI_STATUS_IGNORE,ierr)
|
|
|
|
|
call MPI_File_close(fileUnit,ierr)
|
2019-06-15 19:18:47 +05:30
|
|
|
|
endif restartRead2
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2019-10-24 02:20:01 +05:30
|
|
|
|
call utilities_updateGamma(C_minMaxAvg)
|
|
|
|
|
call utilities_saveReferenceStiffness
|
2012-11-12 19:44:39 +05:30
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
end subroutine grid_mech_spectral_basic_init
|
|
|
|
|
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!> @brief solution for the basic scheme with internal iterations
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-12 10:36:59 +05:30
|
|
|
|
function grid_mech_spectral_basic_solution(incInfoIn,timeinc,timeinc_old,stress_BC,rotation_BC) result(solution)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! input data for solution
|
2019-03-12 16:06:18 +05:30
|
|
|
|
character(len=*), intent(in) :: &
|
|
|
|
|
incInfoIn
|
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
|
timeinc, & !< time increment of current solution
|
|
|
|
|
timeinc_old !< time increment of last successful increment
|
|
|
|
|
type(tBoundaryCondition), intent(in) :: &
|
|
|
|
|
stress_BC
|
2019-12-03 00:36:58 +05:30
|
|
|
|
type(rotation), intent(in) :: &
|
|
|
|
|
rotation_BC
|
2019-03-12 16:06:18 +05:30
|
|
|
|
type(tSolutionState) :: &
|
|
|
|
|
solution
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-08-08 14:43:29 +05:30
|
|
|
|
! PETSc Data
|
2019-03-12 16:06:18 +05:30
|
|
|
|
PetscErrorCode :: ierr
|
|
|
|
|
SNESConvergedReason :: reason
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
incInfo = incInfoIn
|
2012-12-14 20:48:04 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! update stiffness (and gamma operator)
|
2019-04-15 19:23:46 +05:30
|
|
|
|
S = utilities_maskedCompliance(rotation_BC,stress_BC%maskLogical,C_volAvg)
|
2019-10-25 04:12:59 +05:30
|
|
|
|
if(num%update_gamma) call utilities_updateGamma(C_minMaxAvg)
|
2019-10-24 02:20:01 +05:30
|
|
|
|
|
2013-01-10 19:03:43 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-25 23:23:12 +05:30
|
|
|
|
! set module wide available data
|
2019-03-12 16:06:18 +05:30
|
|
|
|
params%stress_mask = stress_BC%maskFloat
|
|
|
|
|
params%stress_BC = stress_BC%values
|
|
|
|
|
params%rotation_BC = rotation_BC
|
|
|
|
|
params%timeinc = timeinc
|
|
|
|
|
params%timeincOld = timeinc_old
|
2012-11-12 19:44:39 +05:30
|
|
|
|
|
2015-12-02 04:06:19 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-25 23:23:12 +05:30
|
|
|
|
! solve BVP
|
2019-03-12 16:06:18 +05:30
|
|
|
|
call SNESsolve(snes,PETSC_NULL_VEC,solution_vec,ierr); CHKERRQ(ierr)
|
2015-12-02 04:06:19 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! check convergence
|
2019-03-12 16:06:18 +05:30
|
|
|
|
call SNESGetConvergedReason(snes,reason,ierr); CHKERRQ(ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
solution%converged = reason > 0
|
|
|
|
|
solution%iterationsNeeded = totalIter
|
|
|
|
|
solution%termIll = terminallyIll
|
|
|
|
|
terminallyIll = .false.
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
end function grid_mech_spectral_basic_solution
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
2013-01-10 03:49:32 +05:30
|
|
|
|
|
2013-12-18 15:05:05 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief forwarding routine
|
2018-02-16 20:06:18 +05:30
|
|
|
|
!> @details find new boundary conditions and best F estimate for end of current timestep
|
|
|
|
|
!> possibly writing restart information, triggering of state increment in DAMASK, and updating of IPcoordinates
|
2013-12-18 15:05:05 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-10-30 03:54:12 +05:30
|
|
|
|
subroutine grid_mech_spectral_basic_forward(cutBack,guess,timeinc,timeinc_old,loadCaseTime,&
|
|
|
|
|
deformation_BC,stress_BC,rotation_BC)
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2018-07-04 00:56:11 +05:30
|
|
|
|
logical, intent(in) :: &
|
2019-10-30 03:54:12 +05:30
|
|
|
|
cutBack, &
|
2018-02-16 20:06:18 +05:30
|
|
|
|
guess
|
2018-07-04 00:56:11 +05:30
|
|
|
|
real(pReal), intent(in) :: &
|
2018-02-16 20:06:18 +05:30
|
|
|
|
timeinc_old, &
|
|
|
|
|
timeinc, &
|
|
|
|
|
loadCaseTime !< remaining time of current load case
|
|
|
|
|
type(tBoundaryCondition), intent(in) :: &
|
|
|
|
|
stress_BC, &
|
|
|
|
|
deformation_BC
|
2019-12-03 00:36:58 +05:30
|
|
|
|
type(rotation), intent(in) :: &
|
2018-02-16 20:06:18 +05:30
|
|
|
|
rotation_BC
|
2019-03-12 10:36:59 +05:30
|
|
|
|
PetscErrorCode :: ierr
|
2018-02-16 20:06:18 +05:30
|
|
|
|
PetscScalar, dimension(:,:,:,:), pointer :: F
|
|
|
|
|
|
|
|
|
|
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2018-02-16 20:06:18 +05:30
|
|
|
|
if (cutBack) then
|
2019-10-24 09:46:42 +05:30
|
|
|
|
C_volAvg = C_volAvgLastInc
|
|
|
|
|
C_minMaxAvg = C_minMaxAvgLastInc
|
2018-02-16 20:06:18 +05:30
|
|
|
|
else
|
|
|
|
|
C_volAvgLastInc = C_volAvg
|
|
|
|
|
C_minMaxAvgLastInc = C_minMaxAvg
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2018-02-25 17:18:58 +05:30
|
|
|
|
F_aimDot = merge(stress_BC%maskFloat*(F_aim-F_aim_lastInc)/timeinc_old, 0.0_pReal, guess)
|
2018-02-16 20:06:18 +05:30
|
|
|
|
F_aim_lastInc = F_aim
|
2018-02-25 17:18:58 +05:30
|
|
|
|
|
2020-01-03 18:32:14 +05:30
|
|
|
|
!-----------------------------------------------------------------------------------------------
|
2018-02-16 20:06:18 +05:30
|
|
|
|
! calculate rate for aim
|
2019-10-24 02:36:47 +05:30
|
|
|
|
if (deformation_BC%myType=='l') then ! calculate F_aimDot from given L and current F
|
2018-02-16 20:06:18 +05:30
|
|
|
|
F_aimDot = &
|
2019-04-03 11:52:04 +05:30
|
|
|
|
F_aimDot + deformation_BC%maskFloat * matmul(deformation_BC%values, F_aim_lastInc)
|
2019-10-24 02:36:47 +05:30
|
|
|
|
elseif(deformation_BC%myType=='fdot') then ! F_aimDot is prescribed
|
2018-02-16 20:06:18 +05:30
|
|
|
|
F_aimDot = &
|
|
|
|
|
F_aimDot + deformation_BC%maskFloat * deformation_BC%values
|
2019-10-24 02:36:47 +05:30
|
|
|
|
elseif (deformation_BC%myType=='f') then ! aim at end of load case is prescribed
|
2018-02-16 20:06:18 +05:30
|
|
|
|
F_aimDot = &
|
|
|
|
|
F_aimDot + deformation_BC%maskFloat * (deformation_BC%values - F_aim_lastInc)/loadCaseTime
|
|
|
|
|
endif
|
|
|
|
|
|
2019-12-03 00:36:58 +05:30
|
|
|
|
Fdot = utilities_calculateRate(guess, &
|
|
|
|
|
F_lastInc,reshape(F,[3,3,grid(1),grid(2),grid3]),timeinc_old, &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
rotation_BC%rotate(F_aimDot,active=.true.))
|
2019-10-30 03:45:02 +05:30
|
|
|
|
F_lastInc = reshape(F,[3,3,grid(1),grid(2),grid3])
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-30 03:45:02 +05:30
|
|
|
|
materialpoint_F0 = reshape(F, [3,3,1,product(grid(1:2))*grid3])
|
2018-02-16 20:06:18 +05:30
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! update average and local deformation gradients
|
|
|
|
|
F_aim = F_aim_lastInc + F_aimDot * timeinc
|
|
|
|
|
F = reshape(Utilities_forwardField(timeinc,F_lastInc,Fdot, & ! estimate of F at end of time+timeinc that matches rotated F_aim on average
|
2020-02-25 23:23:12 +05:30
|
|
|
|
rotation_BC%rotate(F_aim,active=.true.)),[9,grid(1),grid(2),grid3])
|
2018-02-16 20:06:18 +05:30
|
|
|
|
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
end subroutine grid_mech_spectral_basic_forward
|
2013-10-23 20:52:12 +05:30
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
|
2019-10-25 10:46:03 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief Age
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-07 19:54:45 +05:30
|
|
|
|
subroutine grid_mech_spectral_basic_updateCoords
|
2019-10-25 10:46:03 +05:30
|
|
|
|
|
|
|
|
|
PetscErrorCode :: ierr
|
|
|
|
|
PetscScalar, dimension(:,:,:,:), pointer :: F
|
|
|
|
|
|
|
|
|
|
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
|
|
|
|
call utilities_updateCoords(F)
|
|
|
|
|
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
|
|
|
|
|
2019-10-28 18:06:36 +05:30
|
|
|
|
end subroutine grid_mech_spectral_basic_updateCoords
|
2019-10-25 10:46:03 +05:30
|
|
|
|
|
|
|
|
|
|
2019-10-24 17:16:36 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief Write current solver and constitutive data for restart to file
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-07 19:54:45 +05:30
|
|
|
|
subroutine grid_mech_spectral_basic_restartWrite
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
|
|
|
|
PetscErrorCode :: ierr
|
2019-12-07 19:54:45 +05:30
|
|
|
|
integer(HID_T) :: fileHandle, groupHandle
|
2019-10-24 17:16:36 +05:30
|
|
|
|
PetscScalar, dimension(:,:,:,:), pointer :: F
|
2019-12-07 19:54:45 +05:30
|
|
|
|
character(len=pStringLen) :: fileName
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
|
|
|
|
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
|
|
|
|
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(6,'(a)') ' writing solver data required for restart to file'; flush(6)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-12-07 19:54:45 +05:30
|
|
|
|
write(fileName,'(a,a,i0,a)') trim(getSolverJobName()),'_',worldrank,'.hdf5'
|
|
|
|
|
fileHandle = HDF5_openFile(fileName,'w')
|
2019-10-28 17:47:05 +05:30
|
|
|
|
groupHandle = HDF5_addGroup(fileHandle,'solver')
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_write(groupHandle,F_aim, 'F_aim')
|
|
|
|
|
call HDF5_write(groupHandle,F_aim_lastInc,'F_aim_lastInc')
|
|
|
|
|
call HDF5_write(groupHandle,F_aimDot, 'F_aimDot')
|
|
|
|
|
call HDF5_write(groupHandle,F, 'F')
|
|
|
|
|
call HDF5_write(groupHandle,F_lastInc, 'F_lastInc')
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_write(groupHandle,C_volAvg, 'C_volAvg')
|
|
|
|
|
call HDF5_write(groupHandle,C_volAvgLastInc,'C_volAvgLastInc')
|
|
|
|
|
call HDF5_write(groupHandle,C_minMaxAvg, 'C_minMaxAvg')
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
2019-10-28 17:47:05 +05:30
|
|
|
|
call HDF5_closeGroup(groupHandle)
|
2019-10-24 17:16:36 +05:30
|
|
|
|
call HDF5_closeFile(fileHandle)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-10-25 02:20:30 +05:30
|
|
|
|
if (num%update_gamma) call utilities_saveReferenceStiffness
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
|
|
|
|
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
2020-01-03 18:32:14 +05:30
|
|
|
|
|
2019-10-24 17:16:36 +05:30
|
|
|
|
end subroutine grid_mech_spectral_basic_restartWrite
|
|
|
|
|
|
|
|
|
|
|
2019-03-23 15:16:56 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief convergence check
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-15 19:23:46 +05:30
|
|
|
|
subroutine converged(snes_local,PETScIter,devNull1,devNull2,devNull3,reason,dummy,ierr)
|
2019-03-23 15:16:56 +05:30
|
|
|
|
|
|
|
|
|
SNES :: snes_local
|
2019-04-15 19:23:46 +05:30
|
|
|
|
PetscInt, intent(in) :: PETScIter
|
|
|
|
|
PetscReal, intent(in) :: &
|
|
|
|
|
devNull1, &
|
|
|
|
|
devNull2, &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
devNull3
|
2019-03-23 15:16:56 +05:30
|
|
|
|
SNESConvergedReason :: reason
|
|
|
|
|
PetscObject :: dummy
|
|
|
|
|
PetscErrorCode :: ierr
|
|
|
|
|
real(pReal) :: &
|
|
|
|
|
divTol, &
|
|
|
|
|
BCTol
|
|
|
|
|
|
|
|
|
|
divTol = max(maxval(abs(P_av))*err_div_tolRel ,err_div_tolAbs)
|
|
|
|
|
BCTol = max(maxval(abs(P_av))*err_stress_tolRel,err_stress_tolAbs)
|
|
|
|
|
|
|
|
|
|
if ((totalIter >= itmin .and. &
|
|
|
|
|
all([ err_div/divTol, &
|
|
|
|
|
err_BC /BCTol ] < 1.0_pReal)) &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
.or. terminallyIll) then
|
2019-03-23 15:16:56 +05:30
|
|
|
|
reason = 1
|
|
|
|
|
elseif (totalIter >= itmax) then
|
|
|
|
|
reason = -1
|
|
|
|
|
else
|
|
|
|
|
reason = 0
|
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! report
|
|
|
|
|
write(6,'(1/,a)') ' ... reporting .............................................................'
|
|
|
|
|
write(6,'(1/,a,f12.2,a,es8.2,a,es9.2,a)') ' error divergence = ', &
|
|
|
|
|
err_div/divTol, ' (',err_div,' / m, tol = ',divTol,')'
|
|
|
|
|
write(6,'(a,f12.2,a,es8.2,a,es9.2,a)') ' error stress BC = ', &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
err_BC/BCTol, ' (',err_BC, ' Pa, tol = ',BCTol,')'
|
2019-03-23 15:16:56 +05:30
|
|
|
|
write(6,'(/,a)') ' ==========================================================================='
|
2020-02-25 23:23:12 +05:30
|
|
|
|
flush(6)
|
|
|
|
|
|
2019-03-23 15:16:56 +05:30
|
|
|
|
end subroutine converged
|
|
|
|
|
|
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-15 19:46:14 +05:30
|
|
|
|
!> @brief forms the residual vector
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
subroutine formResidual(in, F, &
|
|
|
|
|
residuum, dummy, ierr)
|
|
|
|
|
|
|
|
|
|
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: in !< DMDA info (needs to be named "in" for macros like XRANGE to work)
|
|
|
|
|
PetscScalar, dimension(3,3,XG_RANGE,YG_RANGE,ZG_RANGE), &
|
|
|
|
|
intent(in) :: F !< deformation gradient field
|
|
|
|
|
PetscScalar, dimension(3,3,X_RANGE,Y_RANGE,Z_RANGE), &
|
|
|
|
|
intent(out) :: residuum !< residuum field
|
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
|
deltaF_aim
|
|
|
|
|
PetscInt :: &
|
|
|
|
|
PETScIter, &
|
|
|
|
|
nfuncs
|
|
|
|
|
PetscObject :: dummy
|
|
|
|
|
PetscErrorCode :: ierr
|
|
|
|
|
|
|
|
|
|
call SNESGetNumberFunctionEvals(snes,nfuncs,ierr); CHKERRQ(ierr)
|
|
|
|
|
call SNESGetIterationNumber(snes,PETScIter,ierr); CHKERRQ(ierr)
|
|
|
|
|
|
|
|
|
|
if (nfuncs == 0 .and. PETScIter == 0) totalIter = -1 ! new increment
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! begin of new iteration
|
|
|
|
|
newIteration: if (totalIter <= PETScIter) then
|
|
|
|
|
totalIter = totalIter + 1
|
2019-12-11 23:54:29 +05:30
|
|
|
|
write(6,'(1x,a,3(a,i0))') trim(incInfo), ' @ Iteration ', itmin, '≤',totalIter, '≤', itmax
|
2019-03-22 20:32:00 +05:30
|
|
|
|
if (iand(debug_level(debug_spectral),debug_spectralRotation) /= 0) &
|
|
|
|
|
write(6,'(/,a,/,3(3(f12.7,1x)/))',advance='no') &
|
2020-02-25 23:23:12 +05:30
|
|
|
|
' deformation gradient aim (lab) =', transpose(params%rotation_BC%rotate(F_aim,active=.true.))
|
2019-03-22 20:32:00 +05:30
|
|
|
|
write(6,'(/,a,/,3(3(f12.7,1x)/))',advance='no') &
|
|
|
|
|
' deformation gradient aim =', transpose(F_aim)
|
|
|
|
|
flush(6)
|
|
|
|
|
endif newIteration
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! evaluate constitutive response
|
2019-04-15 19:23:46 +05:30
|
|
|
|
call utilities_constitutiveResponse(residuum, & ! "residuum" gets field of first PK stress (to save memory)
|
2019-03-22 20:32:00 +05:30
|
|
|
|
P_av,C_volAvg,C_minMaxAvg, &
|
|
|
|
|
F,params%timeinc,params%rotation_BC)
|
|
|
|
|
call MPI_Allreduce(MPI_IN_PLACE,terminallyIll,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! stress BC handling
|
|
|
|
|
deltaF_aim = math_mul3333xx33(S, P_av - params%stress_BC)
|
|
|
|
|
F_aim = F_aim - deltaF_aim
|
|
|
|
|
err_BC = maxval(abs(params%stress_mask * (P_av - params%stress_BC))) ! mask = 0.0 when no stress bc
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! updated deformation gradient using fix point algorithm of basic scheme
|
|
|
|
|
tensorField_real = 0.0_pReal
|
|
|
|
|
tensorField_real(1:3,1:3,1:grid(1),1:grid(2),1:grid3) = residuum ! store fPK field for subsequent FFT forward transform
|
|
|
|
|
call utilities_FFTtensorForward ! FFT forward of global "tensorField_real"
|
|
|
|
|
err_div = Utilities_divergenceRMS() ! divRMS of tensorField_fourier for later use
|
2020-02-25 23:23:12 +05:30
|
|
|
|
call utilities_fourierGammaConvolution(params%rotation_BC%rotate(deltaF_aim,active=.true.)) ! convolution of Gamma and tensorField_fourier
|
2019-03-22 20:32:00 +05:30
|
|
|
|
call utilities_FFTtensorBackward ! FFT backward of global tensorField_fourier
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! constructing residual
|
|
|
|
|
residuum = tensorField_real(1:3,1:3,1:grid(1),1:grid(2),1:grid3) ! Gamma*P gives correction towards div(P) = 0, so needs to be zero, too
|
|
|
|
|
|
|
|
|
|
end subroutine formResidual
|
|
|
|
|
|
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
end module grid_mech_spectral_basic
|