DAMASK_EICMD/src/grid/DAMASK_grid.f90

485 lines
26 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
use grid_mech_spectral_basic
use grid_mech_spectral_polarisation
use grid_mech_FEM
use grid_damage_spectral
use grid_thermal_spectral
use results
implicit none
!--------------------------------------------------------------------------------------------------
! variables related to information from load case and geom file
real(pReal), dimension(9) :: temp_valueVector = 0.0_pReal !< temporarily from loadcase file when reading in tensors (initialize to 0.0)
logical, dimension(9) :: temp_maskVector = .false. !< temporarily from loadcase file when reading in tensors
2020-06-18 02:30:03 +05:30
integer :: &
N_t = 0, & !< # of time indicators found in load case file
N_n = 0, & !< # of increment specifiers found in load case file
N_def = 0 !< # of rate of deformation specifiers found in load case file
!--------------------------------------------------------------------------------------------------
! 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, &
cutBack = .false.
integer :: &
2020-10-14 14:01:34 +05:30
i, j, m, field, &
2020-06-18 02:30:03 +05:30
errorID = 0, &
cutBackLevel = 0, & !< cut back level \f$ t = \frac{t_{inc}}{2^l} \f$
stepFraction = 0 !< fraction of current time interval
2020-06-18 02:30:03 +05:30
integer :: &
currentLoadcase = 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, &
nActiveFields = 0
character(len=pStringLen) :: &
incInfo, &
loadcase_string
2020-06-29 20:35:11 +05:30
integer :: &
maxCutBack, & !< max number of cut backs
stagItMax !< max number of field level staggered iterations
type(tLoadCase), allocatable, dimension(:) :: loadCases !< array of all load cases
2020-06-18 02:30:03 +05:30
type(tLoadCase) :: newLoadCase
type(tSolutionState), allocatable, dimension(:) :: solres
procedure(grid_mech_spectral_basic_init), pointer :: &
mech_init
procedure(grid_mech_spectral_basic_forward), pointer :: &
mech_forward
procedure(grid_mech_spectral_basic_solution), pointer :: &
mech_solution
procedure(grid_mech_spectral_basic_updateCoords), pointer :: &
mech_updateCoords
procedure(grid_mech_spectral_basic_restartWrite), pointer :: &
mech_restartWrite
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_mech, &
step_discretization, &
step_deformation, &
step_stress
!--------------------------------------------------------------------------------------------------
! init DAMASK (all modules)
2020-06-16 22:17:19 +05:30
2020-06-18 02:30:03 +05:30
call CPFEM_initAll
2020-09-22 16:39:12 +05:30
print'(/,a)', ' <<<+- DAMASK_spectral init -+>>>'; flush(IO_STDOUT)
2020-06-18 02:30:03 +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))
allocate(newLoadCase%ID(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
debug_grid => config_debug%get('grid',defaultVal=emptyList)
2020-06-18 02:30:03 +05:30
select case (trim(num_grid%get_asString('solver', defaultVal = 'Basic')))
case ('Basic')
mech_init => grid_mech_spectral_basic_init
mech_forward => grid_mech_spectral_basic_forward
mech_solution => grid_mech_spectral_basic_solution
mech_updateCoords => grid_mech_spectral_basic_updateCoords
mech_restartWrite => grid_mech_spectral_basic_restartWrite
case ('Polarisation')
2020-06-18 21:13:25 +05:30
if(debug_grid%contains('basic')) &
2020-06-18 02:30:03 +05:30
call IO_warning(42, ext_msg='debug Divergence')
mech_init => grid_mech_spectral_polarisation_init
mech_forward => grid_mech_spectral_polarisation_forward
mech_solution => grid_mech_spectral_polarisation_solution
mech_updateCoords => grid_mech_spectral_polarisation_updateCoords
mech_restartWrite => grid_mech_spectral_polarisation_restartWrite
case ('FEM')
2020-06-18 21:13:25 +05:30
if(debug_grid%contains('basic')) &
2020-06-18 02:30:03 +05:30
call IO_warning(42, ext_msg='debug Divergence')
mech_init => grid_mech_FEM_init
mech_forward => grid_mech_FEM_forward
mech_solution => grid_mech_FEM_solution
mech_updateCoords => grid_mech_FEM_updateCoords
mech_restartWrite => grid_mech_FEM_restartWrite
case default
call IO_error(error_ID = 891, ext_msg = trim(num_grid%get_asString('solver')))
end select
2018-08-31 12:44:16 +05:30
2020-06-18 02:30:03 +05:30
!--------------------------------------------------------------------------------------------------
! reading information from load case file and to sanity checks
2020-10-14 14:01:34 +05:30
config_load => YAML_parse_file(trim(interface_loadFile))
load_steps => config_load%get('step')
allocate(loadCases(load_steps%length)) ! array of load cases
do currentLoadCase = 1, load_steps%length
load_step => load_steps%get(currentLoadCase)
step_mech => load_step%get('mech')
step_discretization => load_step%get('discretization')
do m = 1, step_mech%length
select case (step_mech%getKey(m))
2020-10-16 19:44:13 +05:30
case('L','dot_F','F')
2020-06-18 02:30:03 +05:30
N_def = N_def + 1
end select
enddo
2020-10-14 14:01:34 +05:30
if(step_discretization%contains('t')) N_t = N_t + 1
if(step_discretization%contains('N')) N_n = N_n + 1
if ((N_def /= N_n) .or. (N_n /= N_t) .or. N_n < 1) & ! sanity check
call IO_error(error_ID=837,el=currentLoadCase,ext_msg = trim(interface_loadFile)) ! error message for incomplete loadcase
2020-06-18 02:30:03 +05:30
2020-10-16 03:02:38 +05:30
newLoadCase%stress%myType='P'
2020-06-18 02:30:03 +05:30
field = 1
newLoadCase%ID(field) = FIELD_MECH_ID ! mechanical active by default
2020-06-18 02:30:03 +05:30
thermalActive: if (any(thermal_type == THERMAL_conduction_ID)) then
field = field + 1
newLoadCase%ID(field) = FIELD_THERMAL_ID
endif thermalActive
damageActive: if (any(damage_type == DAMAGE_nonlocal_ID)) then
field = field + 1
newLoadCase%ID(field) = FIELD_DAMAGE_ID
endif damageActive
2020-10-14 14:01:34 +05:30
readMech: do m = 1, step_mech%length
select case (step_mech%getKey(m))
2020-10-16 19:44:13 +05:30
case('dot_F','L','F') ! assign values for the deformation BC matrix
2020-10-14 14:01:34 +05:30
temp_valueVector = 0.0_pReal
2020-10-16 19:44:13 +05:30
if (step_mech%getKey(m) == 'dot_F') then ! in case of dot_F, set type to dot_F
newLoadCase%deformation%myType = 'dot_F'
2020-10-14 14:01:34 +05:30
else if (step_mech%getKey(m) == 'F') then
2020-10-16 03:02:38 +05:30
newLoadCase%deformation%myType = 'F'
2020-06-18 02:30:03 +05:30
else
2020-10-16 03:02:38 +05:30
newLoadCase%deformation%myType = 'L'
2020-06-18 02:30:03 +05:30
endif
2020-10-14 14:01:34 +05:30
step_deformation => step_mech%get(m)
2020-06-18 02:30:03 +05:30
do j = 1, 9
2020-10-16 19:44:13 +05:30
temp_maskVector(j) = step_deformation%get_asString(j) /= 'x' ! true if not a 'x'
2020-10-14 14:01:34 +05:30
if (temp_maskVector(j)) temp_valueVector(j) = step_deformation%get_asFloat(j) ! read value where applicable
2020-06-18 02:30:03 +05:30
enddo
newLoadCase%deformation%mask = transpose(reshape(temp_maskVector,[ 3,3])) ! mask in 3x3 notation
2020-09-24 20:34:06 +05:30
newLoadCase%deformation%values = math_9to33(temp_valueVector) ! values in 3x3 notation
2020-10-14 14:01:34 +05:30
case('P')
2020-06-18 02:30:03 +05:30
temp_valueVector = 0.0_pReal
2020-10-14 14:01:34 +05:30
step_stress => step_mech%get(m)
2020-06-18 02:30:03 +05:30
do j = 1, 9
2020-10-16 19:44:13 +05:30
temp_maskVector(j) = step_stress%get_asString(j) /= 'x' ! true if not a 'x'
2020-10-14 14:01:34 +05:30
if (temp_maskVector(j)) temp_valueVector(j) = step_stress%get_asFloat(j) ! read value where applicable
2020-06-18 02:30:03 +05:30
enddo
newLoadCase%stress%mask = transpose(reshape(temp_maskVector,[ 3,3]))
2020-09-24 20:34:06 +05:30
newLoadCase%stress%values = math_9to33(temp_valueVector)
2020-06-18 02:30:03 +05:30
end select
call newLoadCase%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
!--------------------------------------------------------------------------------------------------------
! Read discretization parameters
!--------------------------------------------------------------------------------------------------------
newLoadCase%time = step_discretization%get_asFloat('t')
newLoadCase%incs = step_discretization%get_asFloat('N')
newLoadCase%logscale = step_discretization%get_asBool ('log_timestep', defaultVal= .false.)
newLoadCase%outputfrequency = step_discretization%get_asInt ('f_out', defaultVal=1)
newLoadCase%restartfrequency = step_discretization%get_asInt ('f_restart', defaultVal=huge(0))
2020-10-14 21:25:57 +05:30
newLoadCase%followFormerTrajectory = .not. (load_step%get_asBool('drop_guessing',defaultVal=.false.) .or. &
merge(.false.,.true.,currentLoadCase > 1)) ! do not continue to predict deformation along former trajectory
2020-06-18 02:30:03 +05:30
reportAndCheck: if (worldrank == 0) then
write (loadcase_string, '(i0)' ) currentLoadCase
2020-09-19 11:50:29 +05:30
print'(/,a,i0)', ' load case: ', currentLoadCase
if (.not. newLoadCase%followFormerTrajectory) &
2020-09-19 11:50:29 +05:30
print*, ' drop guessing along trajectory'
2020-10-16 03:02:38 +05:30
if (newLoadCase%deformation%myType == 'L') then
2020-06-18 02:30:03 +05:30
do j = 1, 3
if (any(newLoadCase%deformation%mask(j,1:3) .eqv. .true.) .and. &
any(newLoadCase%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
2020-09-19 11:50:29 +05:30
print*, ' velocity gradient:'
2020-10-16 03:02:38 +05:30
else if (newLoadCase%deformation%myType == 'F') then
2020-09-19 11:50:29 +05:30
print*, ' deformation gradient at end of load case:'
2020-06-18 02:30:03 +05:30
else
2020-09-19 11:50:29 +05:30
print*, ' deformation gradient rate:'
2020-06-18 02:30:03 +05:30
endif
do i = 1, 3; do j = 1, 3
if(newLoadCase%deformation%mask(i,j)) then
2020-09-22 16:39:12 +05:30
write(IO_STDOUT,'(2x,f12.7)',advance='no') newLoadCase%deformation%values(i,j)
2020-06-18 02:30:03 +05:30
else
2020-09-22 16:39:12 +05:30
write(IO_STDOUT,'(2x,12a)',advance='no') ' * '
2020-06-18 02:30:03 +05:30
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(newLoadCase%stress%mask .eqv. newLoadCase%deformation%mask)) errorID = 831 ! exclusive or masking only
if (any(newLoadCase%stress%mask .and. transpose(newLoadCase%stress%mask) .and. (math_I3<1))) &
errorID = 838 ! no rotation is allowed by stress BC
2020-09-19 11:50:29 +05:30
print*, ' stress / GPa:'
2020-06-18 02:30:03 +05:30
do i = 1, 3; do j = 1, 3
2020-09-25 18:39:24 +05:30
if(newLoadCase%stress%mask(i,j)) then
2020-09-22 16:39:12 +05:30
write(IO_STDOUT,'(2x,f12.7)',advance='no') newLoadCase%stress%values(i,j)*1e-9_pReal
2020-06-18 02:30:03 +05:30
else
2020-09-22 16:39:12 +05:30
write(IO_STDOUT,'(2x,12a)',advance='no') ' * '
2020-06-18 02:30:03 +05:30
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(abs(matmul(newLoadCase%rot%asMatrix(), &
transpose(newLoadCase%rot%asMatrix()))-math_I3) > &
reshape(spread(tol_math_check,1,9),[ 3,3]))) errorID = 846 ! given rotation matrix contains strain
2020-06-18 02:30:03 +05:30
if (any(dNeq(newLoadCase%rot%asMatrix(), math_I3))) &
2020-09-22 16:39:12 +05:30
write(IO_STDOUT,'(2x,a,/,3(3(3x,f12.7,1x)/))',advance='no') 'rotation of loadframe:',&
2020-06-18 02:30:03 +05:30
transpose(newLoadCase%rot%asMatrix())
if (newLoadCase%time < 0.0_pReal) errorID = 834 ! negative time increment
2020-09-19 11:50:29 +05:30
print'(a,f0.3)', ' time: ', newLoadCase%time
if (newLoadCase%incs < 1) errorID = 835 ! non-positive incs count
2020-09-19 11:50:29 +05:30
print'(a,i0)', ' increments: ', newLoadCase%incs
if (newLoadCase%outputfrequency < 1) errorID = 836 ! non-positive result frequency
2020-09-19 11:50:29 +05:30
print'(a,i0)', ' output frequency: ', newLoadCase%outputfrequency
if (newLoadCase%restartfrequency < 1) errorID = 839 ! non-positive restart frequency
2020-06-18 02:30:03 +05:30
if (newLoadCase%restartfrequency < huge(0)) &
2020-09-19 11:50:29 +05:30
print'(a,i0)', ' restart frequency: ', newLoadCase%restartfrequency
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
2020-10-14 14:01:34 +05:30
loadCases(currentLoadCase) = newLoadCase ! load case is ok, append it
2020-06-18 02:30:03 +05:30
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)
call mech_init
case(FIELD_THERMAL_ID)
call grid_thermal_spectral_init
case(FIELD_DAMAGE_ID)
call grid_damage_spectral_init
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
if (debug_grid%contains('basic')) print'(/,a)', ' header of statistics file written out'
2020-09-22 16:39:12 +05:30
flush(IO_STDOUT)
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
writeUndeformed: if (interface_restartInc < 1) then
print'(/,a)', ' ... writing initial configuration to file ........................'
2020-06-18 02:30:03 +05:30
call CPFEM_results(0,0.0_pReal)
endif writeUndeformed
loadCaseLooping: do currentLoadCase = 1, size(loadCases)
time0 = time ! load case start time
guess = loadCases(currentLoadCase)%followFormerTrajectory ! change of load case? homogeneous guess for the first inc
2020-06-18 02:30:03 +05:30
incLooping: do inc = 1, loadCases(currentLoadCase)%incs
totalIncsCounter = totalIncsCounter + 1
!--------------------------------------------------------------------------------------------------
! forwarding time
timeIncOld = timeinc ! last timeinc that brought former inc to an end
2020-10-14 14:01:34 +05:30
if (.not. loadCases(currentLoadCase)%logscale) then ! linear scale
2020-06-18 02:30:03 +05:30
timeinc = loadCases(currentLoadCase)%time/real(loadCases(currentLoadCase)%incs,pReal)
else
if (currentLoadCase == 1) then ! 1st load case of logarithmic scale
2020-09-26 01:05:47 +05:30
timeinc = loadCases(1)%time*(2.0_pReal**real(max(inc-1,1)-loadCases(1)%incs ,pReal)) ! assume 1st inc is equal to 2nd
else ! not-1st load case of logarithmic scale
2020-06-18 02:30:03 +05:30
timeinc = time0 * &
( (1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc ,pReal)/&
real(loadCases(currentLoadCase)%incs ,pReal))&
-(1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc-1 ,pReal)/&
real(loadCases(currentLoadCase)%incs ,pReal)))
endif
endif
timeinc = timeinc * real(subStepFactor,pReal)**real(-cutBackLevel,pReal) ! depending on cut back level, decrease time step
2020-06-18 02:30:03 +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-06-18 02:30:03 +05:30
subStepLooping: do while (stepFraction < subStepFactor**cutBackLevel)
remainingLoadCaseTime = loadCases(currentLoadCase)%time+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(currentLoadCase)%incs,&
'-', stepFraction,'/',subStepFactor**cutBackLevel,&
' of load case ', currentLoadCase,'/',size(loadCases)
write(incInfo,'(4(a,i0))') &
'Increment ',totalIncsCounter,'/',sum(loadCases%incs),&
'-', 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
select case(loadCases(currentLoadCase)%ID(field))
case(FIELD_MECH_ID)
call mech_forward (&
cutBack,guess,timeinc,timeIncOld,remainingLoadCaseTime, &
deformation_BC = loadCases(currentLoadCase)%deformation, &
stress_BC = loadCases(currentLoadCase)%stress, &
rotation_BC = loadCases(currentLoadCase)%rot)
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
select case(loadCases(currentLoadCase)%ID(field))
case(FIELD_MECH_ID)
solres(field) = mech_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
if (.not. solres(field)%converged) exit ! no solution found
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
2020-06-18 02:30:03 +05:30
call mech_updateCoords
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%converged, solres%iterationsNeeded
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) ! quit
2020-06-18 02:30:03 +05:30
endif
enddo subStepLooping
cutBackLevel = max(0, cutBackLevel - 1) ! try half number of subincs next inc
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-06-18 02:30:03 +05:30
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
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
if (mod(inc,loadCases(currentLoadCase)%restartFrequency) == 0) then
call mech_restartWrite
call CPFEM_restartWrite
endif
endif skipping
enddo incLooping
enddo loadCaseLooping
2017-07-27 19:51:02 +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)
call quit(0) ! no complains ;)
2020-06-18 02:30:03 +05:30
end program DAMASK_grid