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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
module grid_mechanical_spectral_basic
|
2018-05-17 15:34:21 +05:30
|
|
|
|
#include <petsc/finclude/petscsnes.h>
|
|
|
|
|
#include <petsc/finclude/petscdmda.h>
|
2021-07-08 18:51:35 +05:30
|
|
|
|
use PETScDMDA
|
|
|
|
|
use PETScSNES
|
2021-07-09 22:18:25 +05:30
|
|
|
|
#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
use MPI_f08
|
|
|
|
|
#endif
|
2019-06-10 00:50:38 +05:30
|
|
|
|
|
|
|
|
|
use prec
|
2020-09-13 13:49:38 +05:30
|
|
|
|
use parallelization
|
2019-06-10 00:50:38 +05:30
|
|
|
|
use DAMASK_interface
|
2020-11-11 16:49:39 +05:30
|
|
|
|
use IO
|
2021-07-08 17:57:04 +05:30
|
|
|
|
use HDF5
|
2019-06-10 00:50:38 +05:30
|
|
|
|
use HDF5_utilities
|
|
|
|
|
use math
|
2020-12-07 21:56:50 +05:30
|
|
|
|
use rotations
|
2019-06-10 00:50:38 +05:30
|
|
|
|
use spectral_utilities
|
|
|
|
|
use config
|
|
|
|
|
use homogenization
|
2020-03-20 11:28:55 +05:30
|
|
|
|
use discretization_grid
|
2019-06-10 00:50:38 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
implicit none
|
|
|
|
|
private
|
2019-03-23 15:16:56 +05:30
|
|
|
|
|
2020-06-29 18:39:13 +05:30
|
|
|
|
type(tSolutionParams) :: params
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2020-09-20 15:19:20 +05:30
|
|
|
|
type :: tNumerics
|
2019-04-15 19:23:46 +05:30
|
|
|
|
logical :: update_gamma !< update gamma operator with current stiffness
|
2020-06-24 20:15:13 +05:30
|
|
|
|
integer :: &
|
|
|
|
|
itmin, & !< minimum number of iterations
|
|
|
|
|
itmax !< maximum number of iterations
|
|
|
|
|
real(pReal) :: &
|
|
|
|
|
eps_div_atol, & !< absolute tolerance for equilibrium
|
|
|
|
|
eps_div_rtol, & !< relative tolerance for equilibrium
|
|
|
|
|
eps_stress_atol, & !< absolute tolerance for fullfillment of stress BC
|
|
|
|
|
eps_stress_rtol !< relative tolerance for fullfillment of stress BC
|
2019-03-26 12:06:55 +05:30
|
|
|
|
end type tNumerics
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2020-06-29 18:39:13 +05:30
|
|
|
|
type(tNumerics) :: num ! numerics parameters. Better name?
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
2020-09-20 15:19:20 +05:30
|
|
|
|
logical :: debugRotation
|
2020-07-02 02:09:44 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! PETSc data
|
2020-06-29 18:39:13 +05:30
|
|
|
|
DM :: da
|
|
|
|
|
SNES :: snes
|
|
|
|
|
Vec :: solution_vec
|
2012-11-12 19:44:39 +05:30
|
|
|
|
|
2012-08-14 22:28:23 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! common pointwise data
|
2020-06-29 18:39:13 +05:30
|
|
|
|
real(pReal), 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.
|
2020-06-29 18:39:13 +05:30
|
|
|
|
real(pReal), dimension(3,3) :: &
|
2019-03-12 16:06:18 +05:30
|
|
|
|
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
|
2020-09-20 16:59:51 +05:30
|
|
|
|
P_av = 0.0_pReal, & !< average 1st Piola--Kirchhoff stress
|
|
|
|
|
P_aim = 0.0_pReal
|
2020-08-09 10:07:50 +05:30
|
|
|
|
character(len=:), allocatable :: incInfo !< time and increment information
|
2020-09-20 15:19:20 +05:30
|
|
|
|
real(pReal), 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
|
|
|
|
|
2020-06-29 18:39:13 +05:30
|
|
|
|
real(pReal) :: &
|
2019-03-12 16:06:18 +05:30
|
|
|
|
err_BC, & !< deviation from stress BC
|
|
|
|
|
err_div !< RMS of div of P
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2020-06-29 18:39:13 +05:30
|
|
|
|
integer :: &
|
2019-03-12 16:06:18 +05:30
|
|
|
|
totalIter = 0 !< total iteration in current increment
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-12 16:06:18 +05:30
|
|
|
|
public :: &
|
2021-02-09 03:51:53 +05:30
|
|
|
|
grid_mechanical_spectral_basic_init, &
|
|
|
|
|
grid_mechanical_spectral_basic_solution, &
|
|
|
|
|
grid_mechanical_spectral_basic_forward, &
|
|
|
|
|
grid_mechanical_spectral_basic_updateCoords, &
|
|
|
|
|
grid_mechanical_spectral_basic_restartWrite
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
subroutine grid_mechanical_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
|
|
|
|
|
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-09-20 15:41:48 +05:30
|
|
|
|
PetscInt, dimension(0:worldsize-1) :: localK
|
2019-12-07 19:54:45 +05:30
|
|
|
|
integer(HID_T) :: fileHandle, groupHandle
|
2021-07-09 22:18:25 +05:30
|
|
|
|
#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
type(MPI_File) :: fileUnit
|
|
|
|
|
#else
|
|
|
|
|
integer :: fileUnit
|
|
|
|
|
#endif
|
2020-09-20 15:41:48 +05:30
|
|
|
|
class (tNode), pointer :: &
|
|
|
|
|
num_grid, &
|
|
|
|
|
debug_grid
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(/,1x,a)', '<<<+- grid_mechanical_spectral_basic init -+>>>'; flush(IO_STDOUT)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(/,1x,a)', 'P. Eisenlohr et al., International Journal of Plasticity 46:37–53, 2013'
|
|
|
|
|
print'( 1x,a)', 'https://doi.org/10.1016/j.ijplas.2012.09.012'//IO_EOL
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'( 1x,a)', 'P. Shanthraj et al., International Journal of Plasticity 66:31–45, 2015'
|
|
|
|
|
print'( 1x,a)', 'https://doi.org/10.1016/j.ijplas.2014.02.006'
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2020-07-02 02:09:44 +05:30
|
|
|
|
!-------------------------------------------------------------------------------------------------
|
|
|
|
|
! debugging options
|
2020-12-15 12:36:50 +05:30
|
|
|
|
debug_grid => config_debug%get('grid',defaultVal=emptyList)
|
2020-07-02 02:31:37 +05:30
|
|
|
|
debugRotation = debug_grid%contains('rotation')
|
2020-12-15 12:36:50 +05:30
|
|
|
|
|
2020-06-18 02:30:03 +05:30
|
|
|
|
!-------------------------------------------------------------------------------------------------
|
2020-06-24 20:15:13 +05:30
|
|
|
|
! read numerical parameters and do sanity checks
|
2020-09-13 14:09:17 +05:30
|
|
|
|
num_grid => config_numerics%get('grid',defaultVal=emptyDict)
|
2020-06-24 20:15:13 +05:30
|
|
|
|
|
2020-12-15 12:36:50 +05:30
|
|
|
|
num%update_gamma = num_grid%get_asBool ('update_gamma', defaultVal=.false.)
|
|
|
|
|
num%eps_div_atol = num_grid%get_asFloat('eps_div_atol', defaultVal=1.0e-4_pReal)
|
|
|
|
|
num%eps_div_rtol = num_grid%get_asFloat('eps_div_rtol', defaultVal=5.0e-4_pReal)
|
|
|
|
|
num%eps_stress_atol = num_grid%get_asFloat('eps_stress_atol',defaultVal=1.0e3_pReal)
|
|
|
|
|
num%eps_stress_rtol = num_grid%get_asFloat('eps_stress_rtol',defaultVal=1.0e-3_pReal)
|
|
|
|
|
num%itmin = num_grid%get_asInt ('itmin',defaultVal=1)
|
|
|
|
|
num%itmax = num_grid%get_asInt ('itmax',defaultVal=250)
|
2020-06-24 20:15:13 +05:30
|
|
|
|
|
|
|
|
|
if (num%eps_div_atol <= 0.0_pReal) call IO_error(301,ext_msg='eps_div_atol')
|
|
|
|
|
if (num%eps_div_rtol < 0.0_pReal) call IO_error(301,ext_msg='eps_div_rtol')
|
|
|
|
|
if (num%eps_stress_atol <= 0.0_pReal) call IO_error(301,ext_msg='eps_stress_atol')
|
|
|
|
|
if (num%eps_stress_rtol < 0.0_pReal) call IO_error(301,ext_msg='eps_stress_rtol')
|
|
|
|
|
if (num%itmax <= 1) call IO_error(301,ext_msg='itmax')
|
|
|
|
|
if (num%itmin > num%itmax .or. num%itmin < 1) call IO_error(301,ext_msg='itmin')
|
2020-12-15 12:36:50 +05:30
|
|
|
|
|
2019-03-12 10:36:59 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! set default and user defined options for PETSc
|
2021-02-09 03:51:53 +05:30
|
|
|
|
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-mechanical_snes_type ngmres',ierr)
|
2019-03-12 16:06:18 +05:30
|
|
|
|
CHKERRQ(ierr)
|
2020-11-11 16:17:23 +05:30
|
|
|
|
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asString('petsc_options',defaultVal=''),ierr)
|
2019-03-12 16:06:18 +05:30
|
|
|
|
CHKERRQ(ierr)
|
2019-03-12 10:36:59 +05:30
|
|
|
|
|
2012-11-12 19:44:39 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! allocate global fields
|
2020-12-15 12:36:50 +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)
|
2021-02-09 03:51:53 +05:30
|
|
|
|
call SNESSetOptionsPrefix(snes,'mechanical_',ierr);CHKERRQ(ierr)
|
2020-12-15 12:36:50 +05:30
|
|
|
|
localK = 0
|
2020-09-20 15:41:48 +05:30
|
|
|
|
localK(worldrank) = grid3
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Allreduce(MPI_IN_PLACE,localK,worldsize,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,ierr)
|
2019-03-12 16:06:18 +05:30
|
|
|
|
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
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(/,1x,a,i0,a)', 'reading restart data of increment ', interface_restartInc, ' from file'
|
2019-03-24 14:21:47 +05:30
|
|
|
|
|
2021-02-22 18:07:21 +05:30
|
|
|
|
fileHandle = HDF5_openFile(getSolverJobName()//'_restart.hdf5','r')
|
2019-10-28 17:47:05 +05:30
|
|
|
|
groupHandle = HDF5_openGroup(fileHandle,'solver')
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(P_aim,groupHandle,'P_aim',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(P_aim,9,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(F_aim,groupHandle,'F_aim',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(F_aim,9,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(F_aim_lastInc,groupHandle,'F_aim_lastInc',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(F_aim_lastInc,9,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(F_aimDot,groupHandle,'F_aimDot',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(F_aimDot,9,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(F,groupHandle,'F')
|
|
|
|
|
call HDF5_read(F_lastInc,groupHandle,'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])
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if restartRead
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-07-17 15:20:21 +05:30
|
|
|
|
homogenization_F0 = reshape(F_lastInc, [3,3,product(grid(1:2))*grid3]) ! set starting condition for homogenization_mechanical_response
|
2020-09-20 15:41:48 +05:30
|
|
|
|
call utilities_updateCoords(reshape(F,shape(F_lastInc)))
|
2020-12-15 12:36:50 +05:30
|
|
|
|
call utilities_constitutiveResponse(P,P_av,C_volAvg,C_minMaxAvg, & ! stress field, stress avg, global average of stiffness and (min+max)/2
|
2019-03-12 16:06:18 +05:30
|
|
|
|
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
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(1x,a,i0,a)', 'reading more restart data of increment ', interface_restartInc, ' from file'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(C_volAvg,groupHandle,'C_volAvg',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(C_volAvg,81,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_read(C_volAvgLastInc,groupHandle,'C_volAvgLastInc',.false.)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Bcast(C_volAvgLastInc,81,MPI_DOUBLE,0,MPI_COMM_WORLD,ierr)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (ierr /=0) error stop 'MPI error'
|
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)
|
|
|
|
|
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_File_open(MPI_COMM_WORLD, trim(getSolverJobName())//'.C_ref', &
|
2019-09-23 04:04:05 +05:30
|
|
|
|
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)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if 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
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end subroutine grid_mechanical_spectral_basic_init
|
2019-03-12 10:36:59 +05:30
|
|
|
|
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
function grid_mechanical_spectral_basic_solution(incInfoIn) 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
|
|
|
|
|
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)
|
2020-09-20 19:54:08 +05:30
|
|
|
|
S = utilities_maskedCompliance(params%rotation_BC,params%stress_mask,C_volAvg)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (num%update_gamma) call utilities_updateGamma(C_minMaxAvg)
|
2019-10-24 02:20:01 +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.
|
2021-07-20 02:26:34 +05:30
|
|
|
|
P_aim = merge(P_av,P_aim,params%stress_mask)
|
2012-08-14 22:28:23 +05:30
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end function grid_mechanical_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
|
2013-12-18 15:05:05 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
subroutine grid_mechanical_spectral_basic_forward(cutBack,guess,Delta_t,Delta_t_old,t_remaining,&
|
2019-10-30 03:54:12 +05:30
|
|
|
|
deformation_BC,stress_BC,rotation_BC)
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2020-12-15 12:36:50 +05:30
|
|
|
|
logical, intent(in) :: &
|
2019-10-30 03:54:12 +05:30
|
|
|
|
cutBack, &
|
2018-02-16 20:06:18 +05:30
|
|
|
|
guess
|
2020-12-15 12:36:50 +05:30
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
|
Delta_t_old, &
|
|
|
|
|
Delta_t, &
|
|
|
|
|
t_remaining !< remaining time of current load case
|
|
|
|
|
type(tBoundaryCondition), intent(in) :: &
|
2018-02-16 20:06:18 +05:30
|
|
|
|
stress_BC, &
|
|
|
|
|
deformation_BC
|
2020-12-15 12:36:50 +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
|
2020-09-20 15:41:48 +05:30
|
|
|
|
PetscScalar, pointer, dimension(:,:,:,:) :: F
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2020-09-14 18:28:44 +05:30
|
|
|
|
|
2018-02-16 20:06:18 +05:30
|
|
|
|
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
|
|
|
|
|
2021-07-20 02:26:34 +05:30
|
|
|
|
F_aimDot = merge(merge(.0_pReal,(F_aim-F_aim_lastInc)/Delta_t_old,stress_BC%mask),.0_pReal,guess) ! estimate deformation rate for prescribed stress components
|
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
|
2020-10-16 03:02:38 +05:30
|
|
|
|
if (deformation_BC%myType=='L') then ! calculate F_aimDot from given L and current F
|
2020-09-20 19:54:08 +05:30
|
|
|
|
F_aimDot = F_aimDot &
|
more flexibility for the L in the load case
Note that mixed boundary conditions for L introduce an ambiguity.
Consider:
L = [[1.0, x, x],
[ 0, 0, 0],
[ 0, 0, 0]]
P = [[x, 0, 0],
[x, x, x],
[x, x, x]]
What we need is F^(n+1)=F_dot^(n+1) x Delta_t, where F_dot^(n+1) is
F_dot^(n+1)_ij = L_ik F^n_kj.
So component F_11 has contributions from L_12 and L_13. We first assume
L_12=L_13=0 and then choose F^(n+1)_12 and F^(n+1)_13 to get
P_12=P_13=0. This implicitly gives a solution for L_12 and L_13, which
is however only one out of infinitely many.
2021-07-20 10:08:29 +05:30
|
|
|
|
+ matmul(merge(.0_pReal,deformation_BC%values,deformation_BC%mask),F_aim_lastInc)
|
2020-12-15 12:36:50 +05:30
|
|
|
|
elseif (deformation_BC%myType=='dot_F') then ! F_aimDot is prescribed
|
|
|
|
|
F_aimDot = F_aimDot &
|
2021-07-20 02:26:34 +05:30
|
|
|
|
+ merge(.0_pReal,deformation_BC%values,deformation_BC%mask)
|
2020-10-16 03:02:38 +05:30
|
|
|
|
elseif (deformation_BC%myType=='F') then ! aim at end of load case is prescribed
|
2020-09-20 19:54:08 +05:30
|
|
|
|
F_aimDot = F_aimDot &
|
2021-07-20 02:26:34 +05:30
|
|
|
|
+ merge(.0_pReal,(deformation_BC%values - F_aim_lastInc)/t_remaining,deformation_BC%mask)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
2019-12-03 00:36:58 +05:30
|
|
|
|
Fdot = utilities_calculateRate(guess, &
|
2020-12-15 12:36:50 +05:30
|
|
|
|
F_lastInc,reshape(F,[3,3,grid(1),grid(2),grid3]),Delta_t_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
|
|
|
|
|
2020-12-16 17:18:45 +05:30
|
|
|
|
homogenization_F0 = reshape(F,[3,3,product(grid(1:2))*grid3])
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! update average and local deformation gradients
|
2020-12-15 12:36:50 +05:30
|
|
|
|
F_aim = F_aim_lastInc + F_aimDot * Delta_t
|
|
|
|
|
if (stress_BC%myType=='P') P_aim = P_aim &
|
2021-07-20 02:26:34 +05:30
|
|
|
|
+ merge(.0_pReal,(stress_BC%values - P_aim)/t_remaining,stress_BC%mask)*Delta_t
|
2020-12-15 12:36:50 +05:30
|
|
|
|
if (stress_BC%myType=='dot_P') P_aim = P_aim &
|
2021-07-20 02:26:34 +05:30
|
|
|
|
+ merge(.0_pReal,stress_BC%values,stress_BC%mask)*Delta_t
|
2020-12-15 12:36:50 +05:30
|
|
|
|
|
|
|
|
|
F = reshape(utilities_forwardField(Delta_t,F_lastInc,Fdot, & ! estimate of F at end of time+Delta_t 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
|
|
|
|
|
2020-12-15 12:36:50 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! set module wide available data
|
|
|
|
|
params%stress_mask = stress_BC%mask
|
|
|
|
|
params%rotation_BC = rotation_BC
|
2021-07-20 02:00:20 +05:30
|
|
|
|
params%Delta_t = Delta_t
|
2020-12-15 12:36:50 +05:30
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end subroutine grid_mechanical_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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-12-15 12:36:50 +05:30
|
|
|
|
!> @brief Update coordinates
|
2019-10-25 10:46:03 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
subroutine grid_mechanical_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)
|
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end subroutine grid_mechanical_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
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-09 03:51:53 +05:30
|
|
|
|
subroutine grid_mechanical_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
|
|
|
|
|
|
|
|
|
|
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
|
|
|
|
|
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(1x,a)', 'writing solver data required for restart to file'; flush(IO_STDOUT)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2021-02-22 18:07:21 +05:30
|
|
|
|
fileHandle = HDF5_openFile(getSolverJobName()//'_restart.hdf5','w')
|
2019-10-28 17:47:05 +05:30
|
|
|
|
groupHandle = HDF5_addGroup(fileHandle,'solver')
|
2021-06-01 20:16:24 +05:30
|
|
|
|
call HDF5_write(F,groupHandle,'F')
|
|
|
|
|
call HDF5_write(F_lastInc,groupHandle,'F_lastInc')
|
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
|
|
|
|
|
2021-08-15 11:34:27 +05:30
|
|
|
|
if (worldrank == 0) then
|
|
|
|
|
fileHandle = HDF5_openFile(getSolverJobName()//'_restart.hdf5','a',.false.)
|
|
|
|
|
groupHandle = HDF5_openGroup(fileHandle,'solver')
|
|
|
|
|
call HDF5_write(P_aim,groupHandle,'P_aim',.false.)
|
|
|
|
|
call HDF5_write(F_aim,groupHandle,'F_aim',.false.)
|
|
|
|
|
call HDF5_write(F_aim_lastInc,groupHandle,'F_aim_lastInc',.false.)
|
|
|
|
|
call HDF5_write(F_aimDot,groupHandle,'F_aimDot',.false.)
|
|
|
|
|
call HDF5_write(C_volAvg,groupHandle,'C_volAvg',.false.)
|
|
|
|
|
call HDF5_write(C_volAvgLastInc,groupHandle,'C_volAvgLastInc',.false.)
|
|
|
|
|
call HDF5_write(C_minMaxAvg,groupHandle,'C_minMaxAvg',.false.)
|
|
|
|
|
call HDF5_closeGroup(groupHandle)
|
|
|
|
|
call HDF5_closeFile(fileHandle)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if
|
2021-08-15 11:34:27 +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
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end subroutine grid_mechanical_spectral_basic_restartWrite
|
2019-10-24 17:16:36 +05:30
|
|
|
|
|
|
|
|
|
|
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, &
|
2020-06-24 20:15:13 +05:30
|
|
|
|
BCTol
|
2020-06-19 06:02:33 +05:30
|
|
|
|
|
2021-07-17 00:02:21 +05:30
|
|
|
|
divTol = max(maxval(abs(P_av))*num%eps_div_rtol, num%eps_div_atol)
|
|
|
|
|
BCTol = max(maxval(abs(P_av))*num%eps_stress_rtol, num%eps_stress_atol)
|
2020-06-24 20:15:13 +05:30
|
|
|
|
|
2021-03-28 03:47:04 +05:30
|
|
|
|
if ((totalIter >= num%itmin .and. all([err_div/divTol, err_BC/BCTol] < 1.0_pReal)) &
|
|
|
|
|
.or. terminallyIll) then
|
2019-03-23 15:16:56 +05:30
|
|
|
|
reason = 1
|
2020-06-24 20:15:13 +05:30
|
|
|
|
elseif (totalIter >= num%itmax) then
|
2019-03-23 15:16:56 +05:30
|
|
|
|
reason = -1
|
|
|
|
|
else
|
|
|
|
|
reason = 0
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if
|
2019-03-23 15:16:56 +05:30
|
|
|
|
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(/,1x,a)', '... reporting .............................................................'
|
|
|
|
|
print'(/,1x,a,f12.2,a,es8.2,a,es9.2,a)', 'error divergence = ', &
|
2019-03-23 15:16:56 +05:30
|
|
|
|
err_div/divTol, ' (',err_div,' / m, tol = ',divTol,')'
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(1x,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,')'
|
2021-11-15 23:05:44 +05:30
|
|
|
|
print'(/,1x,a)', '==========================================================================='
|
2020-09-22 16:39:12 +05:30
|
|
|
|
flush(IO_STDOUT)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
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)
|
2020-12-15 12:36:50 +05:30
|
|
|
|
call SNESGetIterationNumber(snes,PETScIter,ierr); CHKERRQ(ierr)
|
2019-03-22 20:32:00 +05:30
|
|
|
|
|
|
|
|
|
if (nfuncs == 0 .and. PETScIter == 0) totalIter = -1 ! new increment
|
2020-09-20 15:41:48 +05:30
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! begin of new iteration
|
|
|
|
|
newIteration: if (totalIter <= PETScIter) then
|
|
|
|
|
totalIter = totalIter + 1
|
2020-09-19 11:50:29 +05:30
|
|
|
|
print'(1x,a,3(a,i0))', trim(incInfo), ' @ Iteration ', num%itmin, '≤',totalIter, '≤', num%itmax
|
2021-11-15 23:05:44 +05:30
|
|
|
|
if (debugRotation) print'(/,1x,a,/,2(3(f12.7,1x)/),3(f12.7,1x))', &
|
|
|
|
|
'deformation gradient aim (lab) =', transpose(params%rotation_BC%rotate(F_aim,active=.true.))
|
|
|
|
|
print'(/,1x,a,/,2(3(f12.7,1x)/),3(f12.7,1x))', &
|
|
|
|
|
'deformation gradient aim =', transpose(F_aim)
|
2020-09-22 16:39:12 +05:30
|
|
|
|
flush(IO_STDOUT)
|
2021-11-15 23:05:44 +05:30
|
|
|
|
end if newIteration
|
2019-03-22 20:32:00 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! 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, &
|
2021-07-20 02:00:20 +05:30
|
|
|
|
F,params%Delta_t,params%rotation_BC)
|
2021-07-08 18:51:35 +05:30
|
|
|
|
call MPI_Allreduce(MPI_IN_PLACE,terminallyIll,1,MPI_LOGICAL,MPI_LOR,MPI_COMM_WORLD,ierr)
|
2020-02-25 23:23:12 +05:30
|
|
|
|
|
2019-03-22 20:32:00 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! stress BC handling
|
2020-09-20 16:59:51 +05:30
|
|
|
|
deltaF_aim = math_mul3333xx33(S, P_av - P_aim) ! S = 0.0 for no bc
|
|
|
|
|
F_aim = F_aim - deltaF_aim
|
2021-07-20 02:26:34 +05:30
|
|
|
|
err_BC = maxval(abs(merge(.0_pReal,P_av - P_aim,params%stress_mask)))
|
2019-03-22 20:32:00 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! 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"
|
2020-09-20 15:41:48 +05:30
|
|
|
|
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
|
|
|
|
|
|
|
|
|
|
|
2021-02-09 03:51:53 +05:30
|
|
|
|
end module grid_mechanical_spectral_basic
|