DAMASK_EICMD/src/grid/DAMASK_grid.f90

488 lines
25 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @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
!> @brief Driver controlling inner and outer load case looping of the various spectral solvers
!> @details doing cutbacking, forwarding in case of restart, reporting statistics, writing
!> results
!--------------------------------------------------------------------------------------------------
program DAMASK_grid
#include <petsc/finclude/petscsys.h>
2020-06-18 02:30:03 +05:30
use PETScsys
use prec
2020-09-13 13:49:38 +05:30
use parallelization
2020-06-18 02:30:03 +05:30
use DAMASK_interface
use IO
use config
use math
use CPFEM2
use material
use spectral_utilities
2021-02-09 03:51:53 +05:30
use grid_mechanical_spectral_basic
use grid_mechanical_spectral_polarisation
use grid_mechanical_FEM
2020-06-18 02:30:03 +05:30
use grid_damage_spectral
use grid_thermal_spectral
use results
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
implicit none
type :: tLoadCase
type(rotation) :: rot !< rotation of BC
type(tBoundaryCondition) :: stress, & !< stress BC
deformation !< deformation BC (dot_F, F, or L)
real(pReal) :: t, & !< length of increment
r !< ratio of geometric progression
integer :: N, & !< number of increments
f_out, & !< frequency of result writes
f_restart !< frequency of restart writes
2021-02-05 21:29:28 +05:30
logical :: estimate_rate !< follow trajectory of former loadcase
integer(kind(FIELD_UNDEFINED_ID)), allocatable :: ID(:)
end type tLoadCase
!--------------------------------------------------------------------------------------------------
! variables related to information from load case and geom file
real(pReal), dimension(9) :: temp_valueVector !< temporarily from loadcase file when reading in tensors (initialize to 0.0)
logical, dimension(9) :: temp_maskVector !< temporarily from loadcase file when reading in tensors
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
2020-06-18 02:30:03 +05:30
real(pReal), dimension(3,3), parameter :: &
ones = 1.0_pReal, &
zeros = 0.0_pReal
integer, parameter :: &
subStepFactor = 2 !< for each substep, divide the last time increment by 2.0
2020-06-18 02:30:03 +05:30
real(pReal) :: &
time = 0.0_pReal, & !< elapsed time
time0 = 0.0_pReal, & !< begin of interval
timeinc = 1.0_pReal, & !< current time interval
timeIncOld = 0.0_pReal, & !< previous time interval
remainingLoadCaseTime = 0.0_pReal !< remaining time of current load case
2020-06-18 02:30:03 +05:30
logical :: &
guess, & !< guess along former trajectory
2020-06-18 02:30:03 +05:30
stagIterate, &
2020-12-17 17:12:17 +05:30
cutBack = .false.,&
2020-12-17 17:50:18 +05:30
signal
2020-06-18 02:30:03 +05:30
integer :: &
2020-10-14 14:01:34 +05:30
i, j, m, field, &
2020-06-18 02:30:03 +05:30
errorID = 0, &
2020-12-17 17:12:17 +05:30
ierr,&
cutBackLevel = 0, & !< cut back level \f$ t = \frac{t_{inc}}{2^l} \f$
2020-10-22 00:42:55 +05:30
stepFraction = 0, & !< fraction of current time interval
2020-10-21 20:49:15 +05:30
l = 0, & !< current load case
inc, & !< current increment in current load case
totalIncsCounter = 0, & !< total # of increments
statUnit = 0, & !< file unit for statistics output
2020-06-18 02:30:03 +05:30
stagIter, &
2020-10-22 00:42:55 +05:30
nActiveFields = 0, &
maxCutBack, & !< max number of cut backs
stagItMax !< max number of field level staggered iterations
2020-06-18 02:30:03 +05:30
character(len=pStringLen) :: &
incInfo, &
loadcase_string
2020-10-22 00:42:55 +05:30
type(tLoadCase), allocatable, dimension(:) :: loadCases !< array of all load cases
2020-06-18 02:30:03 +05:30
type(tSolutionState), allocatable, dimension(:) :: solres
2021-02-09 03:51:53 +05:30
procedure(grid_mechanical_spectral_basic_init), pointer :: &
mechanical_init
procedure(grid_mechanical_spectral_basic_forward), pointer :: &
mechanical_forward
procedure(grid_mechanical_spectral_basic_solution), pointer :: &
mechanical_solution
procedure(grid_mechanical_spectral_basic_updateCoords), pointer :: &
mechanical_updateCoords
procedure(grid_mechanical_spectral_basic_restartWrite), pointer :: &
mechanical_restartWrite
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
external :: &
quit
class (tNode), pointer :: &
num_grid, &
2020-10-14 21:25:57 +05:30
debug_grid, & ! pointer to grid debug options
2020-10-14 14:01:34 +05:30
config_load, &
load_steps, &
load_step, &
step_bc, &
2020-10-14 14:01:34 +05:30
step_mech, &
step_discretization, &
step_deformation, &
step_stress
!--------------------------------------------------------------------------------------------------
! init DAMASK (all modules)
2020-10-21 20:49:15 +05:30
call CPFEM_initAll
2020-09-22 16:39:12 +05:30
print'(/,a)', ' <<<+- DAMASK_spectral init -+>>>'; flush(IO_STDOUT)
2020-10-21 20:49:15 +05:30
2020-09-19 11:50:29 +05:30
print*, 'Shanthraj et al., Handbook of Mechanics of Materials, 2019'
print*, 'https://doi.org/10.1007/978-981-10-6855-3_80'
2018-11-18 17:11:05 +05:30
!--------------------------------------------------------------------------------------------------
! initialize field solver information
2020-06-18 02:30:03 +05:30
nActiveFields = 1
if (any(thermal_type == THERMAL_conduction_ID )) nActiveFields = nActiveFields + 1
if (any(damage_type == DAMAGE_nonlocal_ID )) nActiveFields = nActiveFields + 1
allocate(solres(nActiveFields))
2020-06-17 20:17:13 +05:30
!-------------------------------------------------------------------------------------------------
! reading field paramters from numerics file and do sanity checks
num_grid => config_numerics%get('grid', defaultVal=emptyDict)
2020-06-29 20:35:11 +05:30
stagItMax = num_grid%get_asInt('maxStaggeredIter',defaultVal=10)
maxCutBack = num_grid%get_asInt('maxCutBack',defaultVal=3)
2020-06-17 20:17:13 +05:30
2020-06-29 20:35:11 +05:30
if (stagItMax < 0) call IO_error(301,ext_msg='maxStaggeredIter')
if (maxCutBack < 0) call IO_error(301,ext_msg='maxCutBack')
2020-06-17 20:17:13 +05:30
!--------------------------------------------------------------------------------------------------
2020-06-18 02:30:03 +05:30
! assign mechanics solver depending on selected type
2020-10-21 20:49:15 +05:30
debug_grid => config_debug%get('grid',defaultVal=emptyList)
2020-10-21 20:49:15 +05:30
select case (trim(num_grid%get_asString('solver', defaultVal = 'Basic')))
2020-06-18 02:30:03 +05:30
case ('Basic')
2021-02-09 03:51:53 +05:30
mechanical_init => grid_mechanical_spectral_basic_init
mechanical_forward => grid_mechanical_spectral_basic_forward
mechanical_solution => grid_mechanical_spectral_basic_solution
mechanical_updateCoords => grid_mechanical_spectral_basic_updateCoords
mechanical_restartWrite => grid_mechanical_spectral_basic_restartWrite
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
case ('Polarisation')
2021-02-09 03:51:53 +05:30
mechanical_init => grid_mechanical_spectral_polarisation_init
mechanical_forward => grid_mechanical_spectral_polarisation_forward
mechanical_solution => grid_mechanical_spectral_polarisation_solution
mechanical_updateCoords => grid_mechanical_spectral_polarisation_updateCoords
mechanical_restartWrite => grid_mechanical_spectral_polarisation_restartWrite
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
case ('FEM')
2021-02-09 03:51:53 +05:30
mechanical_init => grid_mechanical_FEM_init
mechanical_forward => grid_mechanical_FEM_forward
mechanical_solution => grid_mechanical_FEM_solution
mechanical_updateCoords => grid_mechanical_FEM_updateCoords
mechanical_restartWrite => grid_mechanical_FEM_restartWrite
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
case default
call IO_error(error_ID = 891, ext_msg = trim(num_grid%get_asString('solver')))
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
end select
2020-10-21 20:49:15 +05:30
2018-08-31 12:44:16 +05:30
2020-06-18 02:30:03 +05:30
!--------------------------------------------------------------------------------------------------
2020-10-21 20:49:15 +05:30
! reading information from load case file and to sanity checks
config_load => YAML_parse_file(trim(interface_loadFile))
load_steps => config_load%get('loadstep')
2020-10-14 14:01:34 +05:30
allocate(loadCases(load_steps%length)) ! array of load cases
2020-10-21 20:49:15 +05:30
do l = 1, load_steps%length
allocate(loadCases(l)%ID(nActiveFields))
2020-06-18 02:30:03 +05:30
field = 1
2020-10-21 20:49:15 +05:30
loadCases(l)%ID(field) = FIELD_MECH_ID ! mechanical active by default
thermalActive: if (any(thermal_type == THERMAL_conduction_ID)) then
2020-06-18 02:30:03 +05:30
field = field + 1
2020-10-21 20:49:15 +05:30
loadCases(l)%ID(field) = FIELD_THERMAL_ID
2020-06-18 02:30:03 +05:30
endif thermalActive
damageActive: if (any(damage_type == DAMAGE_nonlocal_ID)) then
2020-06-18 02:30:03 +05:30
field = field + 1
2020-10-21 20:49:15 +05:30
loadCases(l)%ID(field) = FIELD_DAMAGE_ID
2020-06-18 02:30:03 +05:30
endif damageActive
2020-10-21 20:49:15 +05:30
load_step => load_steps%get(l)
step_bc => load_step%get('boundary_conditions')
step_mech => step_bc%get('mechanical')
loadCases(l)%stress%myType=''
2020-10-14 14:01:34 +05:30
readMech: do m = 1, step_mech%length
select case (step_mech%getKey(m))
case ('L','dot_F','F') ! assign values for the deformation BC matrix
2020-10-21 20:49:15 +05:30
loadCases(l)%deformation%myType = step_mech%getKey(m)
2020-10-14 14:01:34 +05:30
step_deformation => step_mech%get(m)
temp_valueVector = 0.0_pReal
2020-06-18 02:30:03 +05:30
do j = 1, 9
temp_maskVector(j) = step_deformation%get_asString(j) /= 'x'
if (temp_maskVector(j)) temp_valueVector(j) = step_deformation%get_asFloat(j)
2020-06-18 02:30:03 +05:30
enddo
loadCases(l)%deformation%mask = transpose(reshape(temp_maskVector,[3,3]))
loadCases(l)%deformation%values = math_9to33(temp_valueVector)
case ('dot_P','P')
loadCases(l)%stress%myType = step_mech%getKey(m)
2020-10-14 14:01:34 +05:30
step_stress => step_mech%get(m)
temp_valueVector = 0.0_pReal
2020-06-18 02:30:03 +05:30
do j = 1, 9
temp_maskVector(j) = step_stress%get_asString(j) /= 'x'
if (temp_maskVector(j)) temp_valueVector(j) = step_stress%get_asFloat(j)
2020-06-18 02:30:03 +05:30
enddo
loadCases(l)%stress%mask = transpose(reshape(temp_maskVector,[3,3]))
2020-10-21 20:49:15 +05:30
loadCases(l)%stress%values = math_9to33(temp_valueVector)
2020-06-18 02:30:03 +05:30
end select
call loadCases(l)%rot%fromAxisAngle(step_mech%get_asFloats('R',defaultVal = real([0.0,0.0,1.0,0.0],pReal)),degrees=.true.)
2020-10-14 14:01:34 +05:30
enddo readMech
if (.not. allocated(loadCases(l)%deformation%myType)) call IO_error(error_ID=837,ext_msg = 'L/dot_F/F missing')
2020-10-21 20:49:15 +05:30
step_discretization => load_step%get('discretization')
if(.not. step_discretization%contains('t')) call IO_error(error_ID=837,ext_msg = 't missing')
if(.not. step_discretization%contains('N')) call IO_error(error_ID=837,ext_msg = 'N missing')
loadCases(l)%t = step_discretization%get_asFloat('t')
loadCases(l)%N = step_discretization%get_asInt ('N')
loadCases(l)%r = step_discretization%get_asFloat('r', defaultVal= 1.0_pReal)
loadCases(l)%f_restart = step_discretization%get_asInt ('f_restart', defaultVal=huge(0))
2020-10-21 20:49:15 +05:30
loadCases(l)%f_out = load_step%get_asInt('f_out', defaultVal=1)
2021-02-05 16:27:18 +05:30
loadCases(l)%estimate_rate = (load_step%get_asBool('estimate_rate',defaultVal=.true.) .and. &
merge(.true.,.false.,l > 1))
2020-10-14 14:01:34 +05:30
2020-06-18 02:30:03 +05:30
reportAndCheck: if (worldrank == 0) then
2020-10-21 20:49:15 +05:30
write (loadcase_string, '(i0)' ) l
print'(/,a,i0)', ' load case: ', l
2021-02-05 16:27:18 +05:30
print*, ' estimate_rate:', loadCases(l)%estimate_rate
2020-10-21 20:49:15 +05:30
if (loadCases(l)%deformation%myType == 'L') then
2020-06-18 02:30:03 +05:30
do j = 1, 3
2020-10-21 20:49:15 +05:30
if (any(loadCases(l)%deformation%mask(j,1:3) .eqv. .true.) .and. &
any(loadCases(l)%deformation%mask(j,1:3) .eqv. .false.)) errorID = 832 ! each row should be either fully or not at all defined
2020-06-18 02:30:03 +05:30
enddo
endif
if (loadCases(l)%deformation%myType == 'F') then
print*, ' F:'
else
print*, ' '//loadCases(l)%deformation%myType//' / 1/s:'
2020-06-18 02:30:03 +05:30
endif
do i = 1, 3; do j = 1, 3
if (loadCases(l)%deformation%mask(i,j)) then
2020-10-21 20:49:15 +05:30
write(IO_STDOUT,'(2x,f12.7)',advance='no') loadCases(l)%deformation%values(i,j)
2020-06-18 02:30:03 +05:30
else
2020-10-21 20:49:15 +05:30
write(IO_STDOUT,'(2x,12a)',advance='no') ' x '
endif
2020-09-22 16:39:12 +05:30
enddo; write(IO_STDOUT,'(/)',advance='no')
2020-06-18 02:30:03 +05:30
enddo
if (any(loadCases(l)%stress%mask .eqv. loadCases(l)%deformation%mask)) errorID = 831
2020-10-21 20:49:15 +05:30
if (any(loadCases(l)%stress%mask .and. transpose(loadCases(l)%stress%mask) .and. (math_I3<1))) &
errorID = 838 ! no rotation is allowed by stress BC
if (loadCases(l)%stress%myType == 'P') print*, ' P / MPa:'
if (loadCases(l)%stress%myType == 'dot_P') print*, ' dot_P / MPa/s:'
if (loadCases(l)%stress%myType /= '') then
do i = 1, 3; do j = 1, 3
if (loadCases(l)%stress%mask(i,j)) then
write(IO_STDOUT,'(2x,f12.4)',advance='no') loadCases(l)%stress%values(i,j)*1e-6_pReal
else
write(IO_STDOUT,'(2x,12a)',advance='no') ' x '
endif
enddo; write(IO_STDOUT,'(/)',advance='no')
enddo
endif
2020-10-21 20:49:15 +05:30
if (any(dNeq(loadCases(l)%rot%asMatrix(), math_I3))) &
write(IO_STDOUT,'(2x,a,/,3(3(3x,f12.7,1x)/))',advance='no') 'R:',&
2020-10-21 20:49:15 +05:30
transpose(loadCases(l)%rot%asMatrix())
if (loadCases(l)%r <= 0.0) errorID = 833
if (loadCases(l)%t < 0.0_pReal) errorID = 834
if (loadCases(l)%N < 1) errorID = 835
if (loadCases(l)%f_out < 1) errorID = 836
if (loadCases(l)%f_restart < 1) errorID = 839
if (dEq(loadCases(l)%r,1.0_pReal,1.e-9_pReal)) then
print'(a)', ' r: 1 (constant step widths)'
else
print'(a,f0.3)', ' r: ', loadCases(l)%r
endif
print'(a,f0.3)', ' t: ', loadCases(l)%t
print'(a,i0)', ' N: ', loadCases(l)%N
print'(a,i0)', ' f_out: ', loadCases(l)%f_out
if (loadCases(l)%f_restart < huge(0)) &
print'(a,i0)', ' f_restart: ', loadCases(l)%f_restart
if (errorID > 0) call IO_error(error_ID = errorID, ext_msg = loadcase_string) ! exit with error message
2020-06-18 02:30:03 +05:30
endif reportAndCheck
enddo
2020-01-21 12:07:04 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-31 13:09:33 +05:30
! doing initialization depending on active solvers
2020-06-18 02:30:03 +05:30
call spectral_Utilities_init
do field = 1, nActiveFields
select case (loadCases(1)%ID(field))
case(FIELD_MECH_ID)
2021-02-09 03:51:53 +05:30
call mechanical_init
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
case(FIELD_THERMAL_ID)
call grid_thermal_spectral_init
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
case(FIELD_DAMAGE_ID)
call grid_damage_spectral_init
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
end select
enddo
!--------------------------------------------------------------------------------------------------
! write header of output file
2020-06-18 02:30:03 +05:30
if (worldrank == 0) then
writeHeader: if (interface_restartInc < 1) then
open(newunit=statUnit,file=trim(getSolverJobName())//'.sta',form='FORMATTED',status='REPLACE')
write(statUnit,'(a)') 'Increment Time CutbackLevel Converged IterationsNeeded' ! statistics file
2020-06-18 02:30:03 +05:30
else writeHeader
open(newunit=statUnit,file=trim(getSolverJobName())//&
'.sta',form='FORMATTED', position='APPEND', status='OLD')
endif writeHeader
endif
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
writeUndeformed: if (interface_restartInc < 1) then
print'(/,a)', ' ... writing initial configuration to file ........................'
flush(IO_STDOUT)
2020-06-18 02:30:03 +05:30
call CPFEM_results(0,0.0_pReal)
endif writeUndeformed
2020-10-21 20:49:15 +05:30
loadCaseLooping: do l = 1, size(loadCases)
time0 = time ! load case start time
2021-02-05 16:27:18 +05:30
guess = loadCases(l)%estimate_rate ! change of load case? homogeneous guess for the first inc
2020-10-21 20:49:15 +05:30
incLooping: do inc = 1, loadCases(l)%N
2020-06-18 02:30:03 +05:30
totalIncsCounter = totalIncsCounter + 1
!--------------------------------------------------------------------------------------------------
! forwarding time
timeIncOld = timeinc ! last timeinc that brought former inc to an end
if (dEq(loadCases(l)%r,1.0_pReal,1.e-9_pReal)) then ! linear scale
timeinc = loadCases(l)%t/real(loadCases(l)%N,pReal)
2020-06-18 02:30:03 +05:30
else
timeinc = loadCases(l)%t * (loadCases(l)%r**(inc-1)-loadCases(l)%r**inc) &
/ (1.0_pReal-loadCases(l)%r**loadCases(l)%N)
2020-06-18 02:30:03 +05:30
endif
timeinc = timeinc * real(subStepFactor,pReal)**real(-cutBackLevel,pReal) ! depending on cut back level, decrease time step
2020-10-21 20:49:15 +05:30
skipping: if (totalIncsCounter <= interface_restartInc) then ! not yet at restart inc?
time = time + timeinc ! just advance time, skip already performed calculation
guess = .true. ! QUESTION:why forced guessing instead of inheriting loadcase preference
2020-06-18 02:30:03 +05:30
else skipping
stepFraction = 0 ! fraction scaled by stepFactor**cutLevel
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
subStepLooping: do while (stepFraction < subStepFactor**cutBackLevel)
remainingLoadCaseTime = loadCases(l)%t+time0 - time
time = time + timeinc ! forward target time
stepFraction = stepFraction + 1 ! count step
!--------------------------------------------------------------------------------------------------
! report begin of new step
2020-09-19 11:50:29 +05:30
print'(/,a)', ' ###########################################################################'
print'(1x,a,es12.5,6(a,i0))', &
2020-06-18 02:30:03 +05:30
'Time', time, &
's: Increment ', inc,'/',loadCases(l)%N,&
2020-06-18 02:30:03 +05:30
'-', stepFraction,'/',subStepFactor**cutBackLevel,&
2020-10-21 20:49:15 +05:30
' of load case ', l,'/',size(loadCases)
2020-06-18 02:30:03 +05:30
write(incInfo,'(4(a,i0))') &
'Increment ',totalIncsCounter,'/',sum(loadCases%N),&
2020-06-18 02:30:03 +05:30
'-', stepFraction,'/',subStepFactor**cutBackLevel
2020-09-22 16:39:12 +05:30
flush(IO_STDOUT)
!--------------------------------------------------------------------------------------------------
! forward fields
2020-06-18 02:30:03 +05:30
do field = 1, nActiveFields
2020-10-21 20:49:15 +05:30
select case(loadCases(l)%ID(field))
2020-06-18 02:30:03 +05:30
case(FIELD_MECH_ID)
2021-02-09 03:51:53 +05:30
call mechanical_forward (&
2020-06-18 02:30:03 +05:30
cutBack,guess,timeinc,timeIncOld,remainingLoadCaseTime, &
2020-10-21 20:49:15 +05:30
deformation_BC = loadCases(l)%deformation, &
stress_BC = loadCases(l)%stress, &
rotation_BC = loadCases(l)%rot)
2020-06-18 02:30:03 +05:30
case(FIELD_THERMAL_ID); call grid_thermal_spectral_forward(cutBack)
case(FIELD_DAMAGE_ID); call grid_damage_spectral_forward(cutBack)
end select
enddo
if(.not. cutBack) call CPFEM_forward
!--------------------------------------------------------------------------------------------------
! solve fields
2020-06-18 02:30:03 +05:30
stagIter = 0
stagIterate = .true.
do while (stagIterate)
do field = 1, nActiveFields
2020-10-21 20:49:15 +05:30
select case(loadCases(l)%ID(field))
2020-06-18 02:30:03 +05:30
case(FIELD_MECH_ID)
2021-02-09 03:51:53 +05:30
solres(field) = mechanical_solution(incInfo)
2020-06-18 02:30:03 +05:30
case(FIELD_THERMAL_ID)
2020-09-25 18:49:31 +05:30
solres(field) = grid_thermal_spectral_solution(timeinc)
2020-06-18 02:30:03 +05:30
case(FIELD_DAMAGE_ID)
2020-09-25 18:49:31 +05:30
solres(field) = grid_damage_spectral_solution(timeinc)
2020-06-18 02:30:03 +05:30
end select
2020-10-21 20:49:15 +05:30
if (.not. solres(field)%converged) exit ! no solution found
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
enddo
stagIter = stagIter + 1
2020-06-29 20:35:11 +05:30
stagIterate = stagIter < stagItMax &
2020-06-18 02:30:03 +05:30
.and. all(solres(:)%converged) &
.and. .not. all(solres(:)%stagConverged) ! stationary with respect to staggered iteration
2020-06-18 02:30:03 +05:30
enddo
!--------------------------------------------------------------------------------------------------
! check solution for either advance or retry
if ( (all(solres(:)%converged .and. solres(:)%stagConverged)) & ! converged
.and. .not. solres(1)%termIll) then ! and acceptable solution found
2021-02-09 03:51:53 +05:30
call mechanical_updateCoords
2020-06-18 02:30:03 +05:30
timeIncOld = timeinc
cutBack = .false.
guess = .true. ! start guessing after first converged (sub)inc
2020-06-18 02:30:03 +05:30
if (worldrank == 0) then
write(statUnit,*) totalIncsCounter, time, cutBackLevel, &
solres(1)%converged, solres(1)%iterationsNeeded
2020-06-18 02:30:03 +05:30
flush(statUnit)
endif
2020-06-29 20:35:11 +05:30
elseif (cutBackLevel < maxCutBack) then ! further cutbacking tolerated?
2020-06-18 02:30:03 +05:30
cutBack = .true.
stepFraction = (stepFraction - 1) * subStepFactor ! adjust to new denominator
2020-06-18 02:30:03 +05:30
cutBackLevel = cutBackLevel + 1
time = time - timeinc ! rewind time
timeinc = timeinc/real(subStepFactor,pReal) ! cut timestep
2020-09-19 11:50:29 +05:30
print'(/,a)', ' cutting back '
else ! no more options to continue
2020-06-18 02:30:03 +05:30
call IO_warning(850)
if (worldrank == 0) close(statUnit)
call quit(0)
2020-06-18 02:30:03 +05:30
endif
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
enddo subStepLooping
2020-10-21 20:49:15 +05:30
cutBackLevel = max(0, cutBackLevel - 1) ! try half number of subincs next inc
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
if (all(solres(:)%converged)) then
2020-09-19 11:50:29 +05:30
print'(/,a,i0,a)', ' increment ', totalIncsCounter, ' converged'
2020-06-18 02:30:03 +05:30
else
2020-09-19 11:50:29 +05:30
print'(/,a,i0,a)', ' increment ', totalIncsCounter, ' NOT converged'
2020-09-22 16:39:12 +05:30
endif; flush(IO_STDOUT)
2020-10-21 20:49:15 +05:30
call MPI_Allreduce(interface_SIGUSR1,signal,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr)
if (ierr /= 0) error stop 'MPI error'
2020-12-17 17:50:18 +05:30
if (mod(inc,loadCases(l)%f_out) == 0 .or. signal) then
2020-09-19 11:50:29 +05:30
print'(1/,a)', ' ... writing results to file ......................................'
2020-09-22 16:39:12 +05:30
flush(IO_STDOUT)
2020-06-18 02:30:03 +05:30
call CPFEM_results(totalIncsCounter,time)
endif
2021-01-18 02:26:19 +05:30
if(signal) call interface_setSIGUSR1(.false.)
call MPI_Allreduce(interface_SIGUSR2,signal,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr)
if (ierr /= 0) error stop 'MPI error'
2020-12-17 17:50:18 +05:30
if (mod(inc,loadCases(l)%f_restart) == 0 .or. signal) then
2021-02-09 03:51:53 +05:30
call mechanical_restartWrite
2020-06-18 02:30:03 +05:30
call CPFEM_restartWrite
endif
2021-01-18 02:26:19 +05:30
if(signal) call interface_setSIGUSR2(.false.)
call MPI_Allreduce(interface_SIGTERM,signal,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr)
if (ierr /= 0) error stop 'MPI error'
2020-12-17 17:50:18 +05:30
if (signal) exit loadCaseLooping
2020-06-18 02:30:03 +05:30
endif skipping
2020-10-21 20:49:15 +05:30
enddo incLooping
2020-10-21 20:49:15 +05:30
2020-06-18 02:30:03 +05:30
enddo loadCaseLooping
2020-10-21 20:49:15 +05:30
!--------------------------------------------------------------------------------------------------
! report summary of whole calculation
2020-09-19 11:50:29 +05:30
print'(/,a)', ' ###########################################################################'
2020-06-18 02:30:03 +05:30
if (worldrank == 0) close(statUnit)
2020-10-21 20:49:15 +05:30
call quit(0) ! no complains ;)
2020-10-21 20:49:15 +05:30
end program DAMASK_grid