2012-07-30 19:36:22 +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
|
|
|
|
!> @brief Driver controlling inner and outer load case looping of the various spectral solvers
|
2013-02-23 05:54:30 +05:30
|
|
|
!> @details doing cutbacking, forwarding in case of restart, reporting statistics, writing
|
|
|
|
!> results
|
2012-07-30 19:36:22 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2016-01-27 22:18:27 +05:30
|
|
|
program DAMASK_spectral
|
2018-05-17 15:34:21 +05:30
|
|
|
#include <petsc/finclude/petscsys.h>
|
|
|
|
use PETScsys
|
2019-05-15 02:14:38 +05:30
|
|
|
use prec
|
|
|
|
use DAMASK_interface
|
|
|
|
use IO
|
|
|
|
use config
|
|
|
|
use debug
|
|
|
|
use math
|
|
|
|
use CPFEM2
|
|
|
|
use material
|
|
|
|
use spectral_utilities
|
2019-03-12 10:36:59 +05:30
|
|
|
use grid_mech_spectral_basic
|
2019-03-12 23:26:28 +05:30
|
|
|
use grid_mech_spectral_polarisation
|
2019-03-23 11:22:26 +05:30
|
|
|
use grid_mech_FEM
|
2019-03-12 04:21:04 +05:30
|
|
|
use grid_damage_spectral
|
2019-03-12 04:07:06 +05:30
|
|
|
use grid_thermal_spectral
|
2018-11-18 17:11:05 +05:30
|
|
|
use results
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2012-07-19 22:54:56 +05:30
|
|
|
implicit none
|
2015-03-25 21:36:19 +05:30
|
|
|
|
2012-07-20 21:03:13 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! variables related to information from load case and geom file
|
2012-11-09 03:03:58 +05:30
|
|
|
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
|
2019-05-15 02:14:38 +05:30
|
|
|
integer, allocatable, dimension(:) :: chunkPos
|
|
|
|
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
|
2019-12-21 17:07:02 +05:30
|
|
|
character(len=pStringLen) :: &
|
2012-07-20 21:03:13 +05:30
|
|
|
line
|
|
|
|
|
2012-07-19 22:54:56 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! loop variables, convergence etc.
|
2012-12-15 21:51:10 +05:30
|
|
|
real(pReal), dimension(3,3), parameter :: &
|
|
|
|
ones = 1.0_pReal, &
|
2016-08-22 19:57:49 +05:30
|
|
|
zeros = 0.0_pReal
|
2019-05-15 02:14:38 +05:30
|
|
|
integer, parameter :: &
|
|
|
|
subStepFactor = 2 !< for each substep, divide the last time increment by 2.0
|
2012-12-15 21:51:10 +05:30
|
|
|
real(pReal) :: &
|
|
|
|
time = 0.0_pReal, & !< elapsed time
|
|
|
|
time0 = 0.0_pReal, & !< begin of interval
|
|
|
|
timeinc = 1.0_pReal, & !< current time interval
|
2013-05-13 15:14:23 +05:30
|
|
|
timeIncOld = 0.0_pReal, & !< previous time interval
|
|
|
|
remainingLoadCaseTime = 0.0_pReal !< remaining time of current load case
|
2012-12-15 21:51:10 +05:30
|
|
|
logical :: &
|
2015-08-21 23:21:05 +05:30
|
|
|
guess, & !< guess along former trajectory
|
2019-10-30 03:54:12 +05:30
|
|
|
stagIterate, &
|
|
|
|
cutBack = .false.
|
2019-05-15 02:14:38 +05:30
|
|
|
integer :: &
|
2015-06-03 23:00:31 +05:30
|
|
|
i, j, k, l, field, &
|
2019-05-15 02:14:38 +05:30
|
|
|
errorID = 0, &
|
|
|
|
cutBackLevel = 0, & !< cut back level \f$ t = \frac{t_{inc}}{2^l} \f$
|
|
|
|
stepFraction = 0 !< fraction of current time interval
|
|
|
|
integer :: &
|
|
|
|
currentLoadcase = 0, & !< current load case
|
2012-12-15 21:51:10 +05:30
|
|
|
inc, & !< current increment in current load case
|
2019-05-15 02:14:38 +05:30
|
|
|
totalIncsCounter = 0, & !< total # of increments
|
|
|
|
fileUnit = 0, & !< file unit for reading load case and writing results
|
2018-08-31 13:44:33 +05:30
|
|
|
myStat, &
|
2019-05-15 02:14:38 +05:30
|
|
|
statUnit = 0, & !< file unit for statistics output
|
2019-10-24 02:08:17 +05:30
|
|
|
stagIter, &
|
|
|
|
nActiveFields = 0
|
2012-07-23 15:42:31 +05:30
|
|
|
character(len=6) :: loadcase_string
|
2018-07-10 11:54:45 +05:30
|
|
|
character(len=1024) :: &
|
2018-09-28 11:19:52 +05:30
|
|
|
incInfo
|
2012-12-15 21:51:10 +05:30
|
|
|
type(tLoadCase), allocatable, dimension(:) :: loadCases !< array of all load cases
|
2018-08-31 12:41:09 +05:30
|
|
|
type(tLoadCase) :: newLoadCase
|
2015-06-03 23:00:31 +05:30
|
|
|
type(tSolutionState), allocatable, dimension(:) :: solres
|
2019-03-12 10:36:59 +05:30
|
|
|
procedure(grid_mech_spectral_basic_init), pointer :: &
|
2018-08-31 11:50:23 +05:30
|
|
|
mech_init
|
2019-03-12 10:36:59 +05:30
|
|
|
procedure(grid_mech_spectral_basic_forward), pointer :: &
|
2018-08-31 11:50:23 +05:30
|
|
|
mech_forward
|
2019-03-12 10:36:59 +05:30
|
|
|
procedure(grid_mech_spectral_basic_solution), pointer :: &
|
2018-08-31 11:50:23 +05:30
|
|
|
mech_solution
|
2019-10-28 18:06:36 +05:30
|
|
|
procedure(grid_mech_spectral_basic_updateCoords), pointer :: &
|
|
|
|
mech_updateCoords
|
2019-10-25 04:12:59 +05:30
|
|
|
procedure(grid_mech_spectral_basic_restartWrite), pointer :: &
|
|
|
|
mech_restartWrite
|
2018-05-17 15:34:21 +05:30
|
|
|
|
2015-08-21 23:21:05 +05:30
|
|
|
external :: &
|
2018-05-17 15:34:21 +05:30
|
|
|
quit
|
|
|
|
|
2012-08-03 14:55:48 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! init DAMASK (all modules)
|
2018-09-24 01:01:30 +05:30
|
|
|
call CPFEM_initAll
|
2019-12-07 19:54:45 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- DAMASK_spectral init -+>>>'; flush(6)
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2019-03-12 04:37:44 +05:30
|
|
|
write(6,'(/,a)') ' Shanthraj et al., Handbook of Mechanics of Materials, 2019'
|
|
|
|
write(6,'(a)') ' https://doi.org/10.1007/978-981-10-6855-3_80'
|
2018-11-18 17:11:05 +05:30
|
|
|
|
2015-06-03 23:00:31 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! initialize field solver information
|
|
|
|
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))
|
2018-08-31 13:09:33 +05:30
|
|
|
allocate(newLoadCase%ID(nActiveFields))
|
2012-10-24 17:01:40 +05:30
|
|
|
|
2018-08-31 11:50:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! assign mechanics solver depending on selected type
|
2019-04-04 03:13:00 +05:30
|
|
|
select case (trim(config_numerics%getString('spectral_solver',defaultVal='basic')))
|
2019-03-23 11:22:26 +05:30
|
|
|
case ('basic')
|
2019-10-24 17:16:36 +05:30
|
|
|
mech_init => grid_mech_spectral_basic_init
|
|
|
|
mech_forward => grid_mech_spectral_basic_forward
|
|
|
|
mech_solution => grid_mech_spectral_basic_solution
|
2019-10-28 18:06:36 +05:30
|
|
|
mech_updateCoords => grid_mech_spectral_basic_updateCoords
|
2019-10-25 04:12:59 +05:30
|
|
|
mech_restartWrite => grid_mech_spectral_basic_restartWrite
|
2018-08-31 11:50:23 +05:30
|
|
|
|
2019-03-23 11:22:26 +05:30
|
|
|
case ('polarisation')
|
2018-08-31 11:50:23 +05:30
|
|
|
if(iand(debug_level(debug_spectral),debug_levelBasic)/= 0) &
|
2019-05-15 02:14:38 +05:30
|
|
|
call IO_warning(42, ext_msg='debug Divergence')
|
2019-10-24 17:16:36 +05:30
|
|
|
mech_init => grid_mech_spectral_polarisation_init
|
|
|
|
mech_forward => grid_mech_spectral_polarisation_forward
|
|
|
|
mech_solution => grid_mech_spectral_polarisation_solution
|
2019-10-28 18:06:36 +05:30
|
|
|
mech_updateCoords => grid_mech_spectral_polarisation_updateCoords
|
2019-10-25 04:12:59 +05:30
|
|
|
mech_restartWrite => grid_mech_spectral_polarisation_restartWrite
|
2019-10-25 11:25:23 +05:30
|
|
|
|
2019-03-23 11:22:26 +05:30
|
|
|
case ('fem')
|
|
|
|
if(iand(debug_level(debug_spectral),debug_levelBasic)/= 0) &
|
2019-05-15 02:14:38 +05:30
|
|
|
call IO_warning(42, ext_msg='debug Divergence')
|
2019-10-24 17:16:36 +05:30
|
|
|
mech_init => grid_mech_FEM_init
|
|
|
|
mech_forward => grid_mech_FEM_forward
|
|
|
|
mech_solution => grid_mech_FEM_solution
|
2019-10-28 18:06:36 +05:30
|
|
|
mech_updateCoords => grid_mech_FEM_updateCoords
|
2019-10-25 04:12:59 +05:30
|
|
|
mech_restartWrite => grid_mech_FEM_restartWrite
|
2018-08-31 11:50:23 +05:30
|
|
|
|
|
|
|
case default
|
2019-05-15 02:14:38 +05:30
|
|
|
call IO_error(error_ID = 891, ext_msg = config_numerics%getString('spectral_solver'))
|
2018-08-31 11:50:23 +05:30
|
|
|
|
|
|
|
end select
|
|
|
|
|
2012-07-19 22:54:56 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-31 13:09:33 +05:30
|
|
|
! reading information from load case file and to sanity checks
|
2018-08-31 13:49:37 +05:30
|
|
|
allocate (loadCases(0)) ! array of load cases
|
2018-08-31 12:28:13 +05:30
|
|
|
open(newunit=fileunit,iostat=myStat,file=trim(loadCaseFile),action='read')
|
2019-05-15 02:14:38 +05:30
|
|
|
if (myStat /= 0) call IO_error(100,el=myStat,ext_msg=trim(loadCaseFile))
|
2012-07-19 22:54:56 +05:30
|
|
|
do
|
2018-08-31 12:28:13 +05:30
|
|
|
read(fileUnit, '(A)', iostat=myStat) line
|
2019-05-15 02:14:38 +05:30
|
|
|
if ( myStat /= 0) exit
|
2012-07-19 22:54:56 +05:30
|
|
|
if (IO_isBlank(line)) cycle ! skip empty lines
|
2018-08-31 13:49:37 +05:30
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
currentLoadCase = currentLoadCase + 1
|
2018-08-31 13:49:37 +05:30
|
|
|
|
2015-08-28 13:08:48 +05:30
|
|
|
chunkPos = IO_stringPos(line)
|
2019-05-15 02:14:38 +05:30
|
|
|
do i = 1, chunkPos(1) ! reading compulsory parameters for loadcase
|
2015-08-28 13:08:48 +05:30
|
|
|
select case (IO_lc(IO_stringValue(line,chunkPos,i)))
|
2019-10-24 10:07:28 +05:30
|
|
|
case('l','fdot','dotf','f')
|
2019-05-15 02:14:38 +05:30
|
|
|
N_def = N_def + 1
|
2013-02-23 05:54:30 +05:30
|
|
|
case('t','time','delta')
|
2019-05-15 02:14:38 +05:30
|
|
|
N_t = N_t + 1
|
2019-10-24 10:07:28 +05:30
|
|
|
case('n','incs','increments','logincs','logincrements')
|
2019-05-15 02:14:38 +05:30
|
|
|
N_n = N_n + 1
|
2013-02-23 05:54:30 +05:30
|
|
|
end select
|
2018-08-31 13:49:37 +05:30
|
|
|
enddo
|
2019-05-15 02:14:38 +05:30
|
|
|
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(loadCaseFile)) ! error message for incomplete loadcase
|
2018-08-31 13:09:33 +05:30
|
|
|
|
|
|
|
newLoadCase%stress%myType='stress'
|
2015-06-03 23:00:31 +05:30
|
|
|
field = 1
|
2018-08-31 13:49:37 +05:30
|
|
|
newLoadCase%ID(field) = FIELD_MECH_ID ! mechanical active by default
|
2016-07-29 20:04:51 +05:30
|
|
|
thermalActive: if (any(thermal_type == THERMAL_conduction_ID)) then
|
2015-06-03 23:00:31 +05:30
|
|
|
field = field + 1
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%ID(field) = FIELD_THERMAL_ID
|
2016-07-29 20:04:51 +05:30
|
|
|
endif thermalActive
|
|
|
|
damageActive: if (any(damage_type == DAMAGE_nonlocal_ID)) then
|
2015-06-03 23:00:31 +05:30
|
|
|
field = field + 1
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%ID(field) = FIELD_DAMAGE_ID
|
2016-07-29 20:04:51 +05:30
|
|
|
endif damageActive
|
2012-07-19 22:54:56 +05:30
|
|
|
|
2019-12-02 21:07:22 +05:30
|
|
|
call newLoadCase%rot%fromEulers(real([0.0,0.0,0.0],pReal))
|
2019-05-15 02:14:38 +05:30
|
|
|
readIn: do i = 1, chunkPos(1)
|
2015-08-28 13:08:48 +05:30
|
|
|
select case (IO_lc(IO_stringValue(line,chunkPos,i)))
|
2019-10-24 10:07:28 +05:30
|
|
|
case('fdot','dotf','l','f') ! assign values for the deformation BC matrix
|
2012-11-09 02:05:31 +05:30
|
|
|
temp_valueVector = 0.0_pReal
|
2015-08-28 13:08:48 +05:30
|
|
|
if (IO_lc(IO_stringValue(line,chunkPos,i)) == 'fdot'.or. & ! in case of Fdot, set type to fdot
|
|
|
|
IO_lc(IO_stringValue(line,chunkPos,i)) == 'dotf') then
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%deformation%myType = 'fdot'
|
2015-08-28 13:08:48 +05:30
|
|
|
else if (IO_lc(IO_stringValue(line,chunkPos,i)) == 'f') then
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%deformation%myType = 'f'
|
2012-11-29 18:56:17 +05:30
|
|
|
else
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%deformation%myType = 'l'
|
2012-08-03 14:55:48 +05:30
|
|
|
endif
|
2019-05-15 02:14:38 +05:30
|
|
|
do j = 1, 9
|
2015-08-28 13:08:48 +05:30
|
|
|
temp_maskVector(j) = IO_stringValue(line,chunkPos,i+j) /= '*' ! true if not a *
|
|
|
|
if (temp_maskVector(j)) temp_valueVector(j) = IO_floatValue(line,chunkPos,i+j) ! read value where applicable
|
2012-07-19 22:54:56 +05:30
|
|
|
enddo
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%deformation%maskLogical = transpose(reshape(temp_maskVector,[ 3,3])) ! logical mask in 3x3 notation
|
|
|
|
newLoadCase%deformation%maskFloat = merge(ones,zeros,newLoadCase%deformation%maskLogical)! float (1.0/0.0) mask in 3x3 notation
|
2019-03-08 02:50:00 +05:30
|
|
|
newLoadCase%deformation%values = math_9to33(temp_valueVector) ! values in 3x3 notation
|
2019-10-24 10:07:28 +05:30
|
|
|
case('p','stress', 's')
|
2012-07-19 22:54:56 +05:30
|
|
|
temp_valueVector = 0.0_pReal
|
2019-05-15 02:14:38 +05:30
|
|
|
do j = 1, 9
|
2015-08-28 13:08:48 +05:30
|
|
|
temp_maskVector(j) = IO_stringValue(line,chunkPos,i+j) /= '*' ! true if not an asterisk
|
|
|
|
if (temp_maskVector(j)) temp_valueVector(j) = IO_floatValue(line,chunkPos,i+j) ! read value where applicable
|
2012-07-19 22:54:56 +05:30
|
|
|
enddo
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%stress%maskLogical = transpose(reshape(temp_maskVector,[ 3,3]))
|
|
|
|
newLoadCase%stress%maskFloat = merge(ones,zeros,newLoadCase%stress%maskLogical)
|
2019-03-08 02:50:00 +05:30
|
|
|
newLoadCase%stress%values = math_9to33(temp_valueVector)
|
2012-07-19 22:54:56 +05:30
|
|
|
case('t','time','delta') ! increment time
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%time = IO_floatValue(line,chunkPos,i+1)
|
2019-10-24 10:07:28 +05:30
|
|
|
case('n','incs','increments') ! number of increments
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%incs = IO_intValue(line,chunkPos,i+1)
|
2019-10-24 10:07:28 +05:30
|
|
|
case('logincs','logincrements') ! number of increments (switch to log time scaling)
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%incs = IO_intValue(line,chunkPos,i+1)
|
|
|
|
newLoadCase%logscale = 1
|
2013-02-23 05:54:30 +05:30
|
|
|
case('freq','frequency','outputfreq') ! frequency of result writings
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%outputfrequency = IO_intValue(line,chunkPos,i+1)
|
2012-07-19 22:54:56 +05:30
|
|
|
case('r','restart','restartwrite') ! frequency of writing restart information
|
2019-06-29 05:39:27 +05:30
|
|
|
newLoadCase%restartfrequency = IO_intValue(line,chunkPos,i+1)
|
2012-07-19 22:54:56 +05:30
|
|
|
case('guessreset','dropguessing')
|
2018-08-31 13:09:33 +05:30
|
|
|
newLoadCase%followFormerTrajectory = .false. ! do not continue to predict deformation along former trajectory
|
|
|
|
case('euler') ! rotation of load case given in euler angles
|
2012-11-09 03:03:58 +05:30
|
|
|
temp_valueVector = 0.0_pReal
|
2019-05-15 02:14:38 +05:30
|
|
|
l = 1 ! assuming values given in degrees
|
|
|
|
k = 1 ! assuming keyword indicating degree/radians present
|
|
|
|
select case (IO_lc(IO_stringValue(line,chunkPos,i+1)))
|
2012-07-19 22:54:56 +05:30
|
|
|
case('deg','degree')
|
2016-08-22 19:57:49 +05:30
|
|
|
case('rad','radian') ! don't convert from degree to radian
|
2019-05-15 02:14:38 +05:30
|
|
|
l = 0
|
2016-08-22 19:57:49 +05:30
|
|
|
case default
|
2019-05-15 02:14:38 +05:30
|
|
|
k = 0
|
2012-07-19 22:54:56 +05:30
|
|
|
end select
|
2019-05-15 02:14:38 +05:30
|
|
|
do j = 1, 3
|
2015-08-28 13:08:48 +05:30
|
|
|
temp_valueVector(j) = IO_floatValue(line,chunkPos,i+k+j)
|
2013-02-06 22:15:34 +05:30
|
|
|
enddo
|
2019-12-02 21:07:22 +05:30
|
|
|
call newLoadCase%rot%fromEulers(temp_valueVector(1:3),degrees=(l==1))
|
2018-08-31 13:09:33 +05:30
|
|
|
case('rotation','rot') ! assign values for the rotation matrix
|
2012-07-19 22:54:56 +05:30
|
|
|
temp_valueVector = 0.0_pReal
|
2019-05-15 02:14:38 +05:30
|
|
|
do j = 1, 9
|
2015-08-28 13:08:48 +05:30
|
|
|
temp_valueVector(j) = IO_floatValue(line,chunkPos,i+j)
|
2013-02-06 22:15:34 +05:30
|
|
|
enddo
|
2019-12-02 21:07:22 +05:30
|
|
|
call newLoadCase%rot%fromMatrix(math_9to33(temp_valueVector))
|
2012-07-19 22:54:56 +05:30
|
|
|
end select
|
2018-08-31 13:44:33 +05:30
|
|
|
enddo readIn
|
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%followFormerTrajectory = merge(.true.,.false.,currentLoadCase > 1) ! by default, guess from previous load case
|
2018-08-31 12:44:16 +05:30
|
|
|
|
2018-08-31 13:09:33 +05:30
|
|
|
reportAndCheck: if (worldrank == 0) then
|
2015-03-25 21:36:19 +05:30
|
|
|
write (loadcase_string, '(i6)' ) currentLoadCase
|
2019-12-13 03:52:37 +05:30
|
|
|
write(6,'(/,1x,a,i0)') 'load case: ', currentLoadCase
|
2018-08-31 13:09:33 +05:30
|
|
|
if (.not. newLoadCase%followFormerTrajectory) write(6,'(2x,a)') 'drop guessing along trajectory'
|
|
|
|
if (newLoadCase%deformation%myType == 'l') then
|
2019-05-15 02:14:38 +05:30
|
|
|
do j = 1, 3
|
2018-08-31 13:09:33 +05:30
|
|
|
if (any(newLoadCase%deformation%maskLogical(j,1:3) .eqv. .true.) .and. &
|
2019-05-15 02:14:38 +05:30
|
|
|
any(newLoadCase%deformation%maskLogical(j,1:3) .eqv. .false.)) errorID = 832 ! each row should be either fully or not at all defined
|
2015-03-25 21:36:19 +05:30
|
|
|
enddo
|
|
|
|
write(6,'(2x,a)') 'velocity gradient:'
|
2018-08-31 13:09:33 +05:30
|
|
|
else if (newLoadCase%deformation%myType == 'f') then
|
2015-03-25 21:36:19 +05:30
|
|
|
write(6,'(2x,a)') 'deformation gradient at end of load case:'
|
|
|
|
else
|
|
|
|
write(6,'(2x,a)') 'deformation gradient rate:'
|
|
|
|
endif
|
2019-05-15 02:14:38 +05:30
|
|
|
do i = 1, 3; do j = 1, 3
|
2018-08-31 13:09:33 +05:30
|
|
|
if(newLoadCase%deformation%maskLogical(i,j)) then
|
|
|
|
write(6,'(2x,f12.7)',advance='no') newLoadCase%deformation%values(i,j)
|
2015-10-26 22:41:36 +05:30
|
|
|
else
|
2016-08-25 21:27:19 +05:30
|
|
|
write(6,'(2x,12a)',advance='no') ' * '
|
2015-10-26 22:41:36 +05:30
|
|
|
endif
|
|
|
|
enddo; write(6,'(/)',advance='no')
|
|
|
|
enddo
|
2018-08-31 13:09:33 +05:30
|
|
|
if (any(newLoadCase%stress%maskLogical .eqv. &
|
2019-05-15 02:14:38 +05:30
|
|
|
newLoadCase%deformation%maskLogical)) errorID = 831 ! exclusive or masking only
|
2019-12-11 04:40:02 +05:30
|
|
|
if (any(newLoadCase%stress%maskLogical .and. transpose(newLoadCase%stress%maskLogical) &
|
|
|
|
.and. (math_I3<1))) errorID = 838 ! no rotation is allowed by stress BC
|
2015-10-26 22:41:36 +05:30
|
|
|
write(6,'(2x,a)') 'stress / GPa:'
|
2019-05-15 02:14:38 +05:30
|
|
|
do i = 1, 3; do j = 1, 3
|
2018-08-31 13:09:33 +05:30
|
|
|
if(newLoadCase%stress%maskLogical(i,j)) then
|
|
|
|
write(6,'(2x,f12.7)',advance='no') newLoadCase%stress%values(i,j)*1e-9_pReal
|
2015-10-26 22:41:36 +05:30
|
|
|
else
|
2016-08-25 21:27:19 +05:30
|
|
|
write(6,'(2x,12a)',advance='no') ' * '
|
2015-10-26 22:41:36 +05:30
|
|
|
endif
|
|
|
|
enddo; write(6,'(/)',advance='no')
|
|
|
|
enddo
|
2019-12-02 21:07:22 +05:30
|
|
|
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
|
|
|
|
if (any(dNeq(newLoadCase%rot%asMatrix(), math_I3))) &
|
2015-03-25 21:36:19 +05:30
|
|
|
write(6,'(2x,a,/,3(3(3x,f12.7,1x)/))',advance='no') 'rotation of loadframe:',&
|
2019-12-02 21:07:22 +05:30
|
|
|
transpose(newLoadCase%rot%asMatrix())
|
2019-05-15 02:14:38 +05:30
|
|
|
if (newLoadCase%time < 0.0_pReal) errorID = 834 ! negative time increment
|
2019-12-13 03:52:37 +05:30
|
|
|
write(6,'(2x,a,f0.3)') 'time: ', newLoadCase%time
|
2019-05-15 02:14:38 +05:30
|
|
|
if (newLoadCase%incs < 1) errorID = 835 ! non-positive incs count
|
2019-12-13 03:52:37 +05:30
|
|
|
write(6,'(2x,a,i0)') 'increments: ', newLoadCase%incs
|
2019-05-15 02:14:38 +05:30
|
|
|
if (newLoadCase%outputfrequency < 1) errorID = 836 ! non-positive result frequency
|
2019-12-13 03:52:37 +05:30
|
|
|
write(6,'(2x,a,i0)') 'output frequency: ', newLoadCase%outputfrequency
|
2019-06-30 03:36:47 +05:30
|
|
|
if (newLoadCase%restartfrequency < 1) errorID = 839 ! non-positive restart frequency
|
2019-06-30 04:59:36 +05:30
|
|
|
if (newLoadCase%restartfrequency < huge(0)) &
|
2019-12-13 03:52:37 +05:30
|
|
|
write(6,'(2x,a,i0)') 'restart frequency: ', newLoadCase%restartfrequency
|
2019-05-15 02:14:38 +05:30
|
|
|
if (errorID > 0) call IO_error(error_ID = errorID, ext_msg = loadcase_string) ! exit with error message
|
2018-08-31 13:09:33 +05:30
|
|
|
endif reportAndCheck
|
|
|
|
loadCases = [loadCases,newLoadCase] ! load case is ok, append it
|
|
|
|
enddo
|
|
|
|
close(fileUnit)
|
2012-07-20 21:03:13 +05:30
|
|
|
|
2012-12-15 21:51:10 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-31 13:09:33 +05:30
|
|
|
! doing initialization depending on active solvers
|
2019-11-24 12:14:17 +05:30
|
|
|
call Utilities_init
|
2015-06-03 23:00:31 +05:30
|
|
|
do field = 1, nActiveFields
|
|
|
|
select case (loadCases(1)%ID(field))
|
|
|
|
case(FIELD_MECH_ID)
|
2018-08-31 11:50:23 +05:30
|
|
|
call mech_init
|
2016-08-25 21:27:19 +05:30
|
|
|
|
2018-08-31 12:41:09 +05:30
|
|
|
case(FIELD_THERMAL_ID)
|
2019-03-12 04:07:06 +05:30
|
|
|
call grid_thermal_spectral_init
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2015-06-03 23:00:31 +05:30
|
|
|
case(FIELD_DAMAGE_ID)
|
2019-03-12 04:21:04 +05:30
|
|
|
call grid_damage_spectral_init
|
2015-06-03 23:00:31 +05:30
|
|
|
|
|
|
|
end select
|
|
|
|
enddo
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2012-07-20 21:03:13 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! write header of output file
|
2015-08-21 23:21:05 +05:30
|
|
|
if (worldrank == 0) then
|
2019-05-15 02:14:38 +05:30
|
|
|
writeHeader: if (interface_restartInc < 1) then
|
2019-09-23 04:26:37 +05:30
|
|
|
open(newunit=statUnit,file=trim(getSolverJobName())//'.sta',form='FORMATTED',status='REPLACE')
|
2015-03-25 21:36:19 +05:30
|
|
|
write(statUnit,'(a)') 'Increment Time CutbackLevel Converged IterationsNeeded' ! statistics file
|
2015-08-21 23:21:05 +05:30
|
|
|
if (iand(debug_level(debug_spectral),debug_levelBasic) /= 0) &
|
2019-12-11 04:40:02 +05:30
|
|
|
write(6,'(/,a)') ' header of statistics file written out'
|
2015-08-21 23:21:05 +05:30
|
|
|
flush(6)
|
2018-08-21 02:44:34 +05:30
|
|
|
else writeHeader
|
2018-07-10 11:54:45 +05:30
|
|
|
open(newunit=statUnit,file=trim(getSolverJobName())//&
|
2015-08-21 23:21:05 +05:30
|
|
|
'.sta',form='FORMATTED', position='APPEND', status='OLD')
|
2018-08-21 02:44:34 +05:30
|
|
|
endif writeHeader
|
2012-07-20 21:03:13 +05:30
|
|
|
endif
|
2013-02-23 05:54:30 +05:30
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
writeUndeformed: if (interface_restartInc < 1) then
|
2018-02-16 20:06:18 +05:30
|
|
|
write(6,'(1/,a)') ' ... writing initial configuration to file ........................'
|
2019-05-15 02:14:38 +05:30
|
|
|
call CPFEM_results(0,0.0_pReal)
|
2018-08-21 02:44:34 +05:30
|
|
|
endif writeUndeformed
|
2018-09-28 11:19:52 +05:30
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
loadCaseLooping: do currentLoadCase = 1, size(loadCases)
|
2018-08-31 13:09:33 +05:30
|
|
|
time0 = time ! load case start time
|
2014-09-04 01:29:47 +05:30
|
|
|
guess = loadCases(currentLoadCase)%followFormerTrajectory ! change of load case? homogeneous guess for the first inc
|
2012-07-19 22:54:56 +05:30
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
incLooping: do inc = 1, loadCases(currentLoadCase)%incs
|
|
|
|
totalIncsCounter = totalIncsCounter + 1
|
2012-07-19 22:54:56 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! forwarding time
|
2018-02-16 20:06:18 +05:30
|
|
|
timeIncOld = timeinc ! last timeinc that brought former inc to an end
|
2019-05-15 02:14:38 +05:30
|
|
|
if (loadCases(currentLoadCase)%logscale == 0) then ! linear scale
|
2018-02-16 20:06:18 +05:30
|
|
|
timeinc = loadCases(currentLoadCase)%time/real(loadCases(currentLoadCase)%incs,pReal)
|
2012-07-19 22:54:56 +05:30
|
|
|
else
|
2019-05-15 02:14:38 +05:30
|
|
|
if (currentLoadCase == 1) then ! 1st load case of logarithmic scale
|
|
|
|
if (inc == 1) then ! 1st inc of 1st load case of logarithmic scale
|
|
|
|
timeinc = loadCases(1)%time*(2.0_pReal**real( 1-loadCases(1)%incs ,pReal)) ! assume 1st inc is equal to 2nd
|
2018-08-31 13:09:33 +05:30
|
|
|
else ! not-1st inc of 1st load case of logarithmic scale
|
2019-05-15 02:14:38 +05:30
|
|
|
timeinc = loadCases(1)%time*(2.0_pReal**real(inc-1-loadCases(1)%incs ,pReal))
|
2012-07-19 22:54:56 +05:30
|
|
|
endif
|
2018-08-31 13:09:33 +05:30
|
|
|
else ! not-1st load case of logarithmic scale
|
2012-12-15 21:51:10 +05:30
|
|
|
timeinc = time0 * &
|
2018-02-17 03:11:07 +05:30
|
|
|
( (1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc ,pReal)/&
|
2012-12-15 21:51:10 +05:30
|
|
|
real(loadCases(currentLoadCase)%incs ,pReal))&
|
2019-05-15 02:14:38 +05:30
|
|
|
-(1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc-1 ,pReal)/&
|
2018-02-16 20:06:18 +05:30
|
|
|
real(loadCases(currentLoadCase)%incs ,pReal)))
|
2012-07-19 22:54:56 +05:30
|
|
|
endif
|
|
|
|
endif
|
2018-02-16 20:06:18 +05:30
|
|
|
timeinc = timeinc * real(subStepFactor,pReal)**real(-cutBackLevel,pReal) ! depending on cut back level, decrease time step
|
2012-07-19 22:54:56 +05:30
|
|
|
|
2019-06-15 19:18:47 +05:30
|
|
|
skipping: if (totalIncsCounter <= interface_restartInc) then ! not yet at restart inc?
|
2018-02-16 20:06:18 +05:30
|
|
|
time = time + timeinc ! just advance time, skip already performed calculation
|
|
|
|
guess = .true. ! QUESTION:why forced guessing instead of inheriting loadcase preference
|
|
|
|
else skipping
|
2019-05-15 02:14:38 +05:30
|
|
|
stepFraction = 0 ! fraction scaled by stepFactor**cutLevel
|
2013-02-23 05:54:30 +05:30
|
|
|
|
2018-02-16 20:06:18 +05:30
|
|
|
subStepLooping: do while (stepFraction < subStepFactor**cutBackLevel)
|
|
|
|
remainingLoadCaseTime = loadCases(currentLoadCase)%time+time0 - time
|
|
|
|
time = time + timeinc ! forward target time
|
2019-05-15 02:14:38 +05:30
|
|
|
stepFraction = stepFraction + 1 ! count step
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2013-12-18 15:05:05 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-02-16 20:06:18 +05:30
|
|
|
! report begin of new step
|
2017-05-24 21:42:36 +05:30
|
|
|
write(6,'(/,a)') ' ###########################################################################'
|
2019-12-07 15:40:48 +05:30
|
|
|
write(6,'(1x,a,es12.5,6(a,i0))') &
|
2017-05-24 21:42:36 +05:30
|
|
|
'Time', time, &
|
2018-02-17 03:11:07 +05:30
|
|
|
's: Increment ', inc,'/',loadCases(currentLoadCase)%incs,&
|
|
|
|
'-', stepFraction,'/',subStepFactor**cutBackLevel,&
|
2017-05-24 21:42:36 +05:30
|
|
|
' of load case ', currentLoadCase,'/',size(loadCases)
|
2019-12-07 15:40:48 +05:30
|
|
|
write(incInfo,'(4(a,i0))') &
|
2018-02-16 20:06:18 +05:30
|
|
|
'Increment ',totalIncsCounter,'/',sum(loadCases%incs),&
|
2018-02-17 03:11:07 +05:30
|
|
|
'-', stepFraction,'/',subStepFactor**cutBackLevel
|
2017-05-24 21:42:36 +05:30
|
|
|
flush(6)
|
2015-06-03 23:00:31 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! forward fields
|
|
|
|
do field = 1, nActiveFields
|
|
|
|
select case(loadCases(currentLoadCase)%ID(field))
|
|
|
|
case(FIELD_MECH_ID)
|
2018-08-31 11:50:23 +05:30
|
|
|
call mech_forward (&
|
2019-10-30 03:54:12 +05:30
|
|
|
cutBack,guess,timeinc,timeIncOld,remainingLoadCaseTime, &
|
2016-10-26 02:24:32 +05:30
|
|
|
deformation_BC = loadCases(currentLoadCase)%deformation, &
|
|
|
|
stress_BC = loadCases(currentLoadCase)%stress, &
|
2019-12-03 00:36:58 +05:30
|
|
|
rotation_BC = loadCases(currentLoadCase)%rot)
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2019-10-30 03:54:12 +05:30
|
|
|
case(FIELD_THERMAL_ID); call grid_thermal_spectral_forward(cutBack)
|
|
|
|
case(FIELD_DAMAGE_ID); call grid_damage_spectral_forward(cutBack)
|
2015-06-03 23:00:31 +05:30
|
|
|
end select
|
2016-08-22 19:57:49 +05:30
|
|
|
enddo
|
2019-10-30 03:45:02 +05:30
|
|
|
if(.not. cutBack) call CPFEM_forward
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2012-12-15 21:51:10 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-06-03 23:00:31 +05:30
|
|
|
! solve fields
|
2019-05-15 02:14:38 +05:30
|
|
|
stagIter = 0
|
2015-06-03 23:00:31 +05:30
|
|
|
stagIterate = .true.
|
|
|
|
do while (stagIterate)
|
|
|
|
do field = 1, nActiveFields
|
|
|
|
select case(loadCases(currentLoadCase)%ID(field))
|
|
|
|
case(FIELD_MECH_ID)
|
2018-08-31 11:50:23 +05:30
|
|
|
solres(field) = mech_solution (&
|
|
|
|
incInfo,timeinc,timeIncOld, &
|
2019-09-23 19:20:25 +05:30
|
|
|
stress_BC = loadCases(currentLoadCase)%stress, &
|
2019-12-03 00:36:58 +05:30
|
|
|
rotation_BC = loadCases(currentLoadCase)%rot)
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2015-06-03 23:00:31 +05:30
|
|
|
case(FIELD_THERMAL_ID)
|
2019-09-23 19:20:25 +05:30
|
|
|
solres(field) = grid_thermal_spectral_solution(timeinc,timeIncOld)
|
2016-08-22 19:57:49 +05:30
|
|
|
|
2015-06-03 23:00:31 +05:30
|
|
|
case(FIELD_DAMAGE_ID)
|
2019-09-23 19:20:25 +05:30
|
|
|
solres(field) = grid_damage_spectral_solution(timeinc,timeIncOld)
|
2015-06-03 23:00:31 +05:30
|
|
|
|
|
|
|
end select
|
2018-02-16 20:06:18 +05:30
|
|
|
|
2016-08-25 21:27:19 +05:30
|
|
|
if (.not. solres(field)%converged) exit ! no solution found
|
2018-02-16 20:06:18 +05:30
|
|
|
|
2015-06-03 23:00:31 +05:30
|
|
|
enddo
|
2019-05-15 02:14:38 +05:30
|
|
|
stagIter = stagIter + 1
|
2018-02-16 20:06:18 +05:30
|
|
|
stagIterate = stagIter < stagItMax &
|
|
|
|
.and. all(solres(:)%converged) &
|
|
|
|
.and. .not. all(solres(:)%stagConverged) ! stationary with respect to staggered iteration
|
2016-08-22 19:57:49 +05:30
|
|
|
enddo
|
2013-02-23 05:54:30 +05:30
|
|
|
|
2012-12-15 21:51:10 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-02-16 20:06:18 +05:30
|
|
|
! check solution for either advance or retry
|
|
|
|
|
2019-10-24 10:02:46 +05:30
|
|
|
if ( (all(solres(:)%converged .and. solres(:)%stagConverged)) & ! converged
|
2018-02-16 20:06:18 +05:30
|
|
|
.and. .not. solres(1)%termIll) then ! and acceptable solution found
|
2019-10-28 18:06:36 +05:30
|
|
|
call mech_updateCoords
|
2018-02-16 20:06:18 +05:30
|
|
|
timeIncOld = timeinc
|
|
|
|
cutBack = .false.
|
2012-11-09 01:02:00 +05:30
|
|
|
guess = .true. ! start guessing after first converged (sub)inc
|
2015-03-25 21:36:19 +05:30
|
|
|
if (worldrank == 0) then
|
|
|
|
write(statUnit,*) totalIncsCounter, time, cutBackLevel, &
|
2018-02-16 20:06:18 +05:30
|
|
|
solres%converged, solres%iterationsNeeded
|
2015-03-25 21:36:19 +05:30
|
|
|
flush(statUnit)
|
2016-08-22 19:57:49 +05:30
|
|
|
endif
|
2018-02-16 20:06:18 +05:30
|
|
|
elseif (cutBackLevel < maxCutBack) then ! further cutbacking tolerated?
|
|
|
|
cutBack = .true.
|
2019-05-15 02:14:38 +05:30
|
|
|
stepFraction = (stepFraction - 1) * subStepFactor ! adjust to new denominator
|
|
|
|
cutBackLevel = cutBackLevel + 1
|
2018-02-16 20:06:18 +05:30
|
|
|
time = time - timeinc ! rewind time
|
|
|
|
timeinc = timeinc/real(subStepFactor,pReal) ! cut timestep
|
|
|
|
write(6,'(/,a)') ' cutting back '
|
|
|
|
else ! no more options to continue
|
2019-05-15 02:14:38 +05:30
|
|
|
call IO_warning(850)
|
2018-02-16 20:06:18 +05:30
|
|
|
close(statUnit)
|
2019-10-25 04:12:59 +05:30
|
|
|
call quit(0) ! quit
|
2016-08-22 19:57:49 +05:30
|
|
|
endif
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
|
|
enddo subStepLooping
|
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
cutBackLevel = max(0, cutBackLevel - 1) ! try half number of subincs next inc
|
2018-02-16 20:06:18 +05:30
|
|
|
|
|
|
|
if (all(solres(:)%converged)) then
|
2019-12-07 19:54:45 +05:30
|
|
|
write(6,'(/,a,i0,a)') ' increment ', totalIncsCounter, ' converged'
|
2012-07-20 21:03:13 +05:30
|
|
|
else
|
2019-12-07 19:54:45 +05:30
|
|
|
write(6,'(/,a,i0,a)') ' increment ', totalIncsCounter, ' NOT converged'
|
2013-02-23 05:54:30 +05:30
|
|
|
endif; flush(6)
|
2018-02-16 20:06:18 +05:30
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
|
2018-02-16 20:06:18 +05:30
|
|
|
write(6,'(1/,a)') ' ... writing results to file ......................................'
|
|
|
|
flush(6)
|
2018-12-17 20:45:16 +05:30
|
|
|
call CPFEM_results(totalIncsCounter,time)
|
2012-07-19 22:54:56 +05:30
|
|
|
endif
|
2019-10-25 12:00:12 +05:30
|
|
|
if (mod(inc,loadCases(currentLoadCase)%restartFrequency) == 0) then
|
|
|
|
call mech_restartWrite
|
|
|
|
call CPFEM_restartWrite
|
|
|
|
endif
|
2018-02-16 20:06:18 +05:30
|
|
|
endif skipping
|
2012-07-20 21:03:13 +05:30
|
|
|
|
2012-07-24 22:37:10 +05:30
|
|
|
enddo incLooping
|
2018-02-16 20:06:18 +05:30
|
|
|
|
2012-07-24 22:37:10 +05:30
|
|
|
enddo loadCaseLooping
|
2017-07-27 19:51:02 +05:30
|
|
|
|
|
|
|
|
2015-03-25 21:36:19 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! report summary of whole calculation
|
2017-05-24 21:42:36 +05:30
|
|
|
write(6,'(/,a)') ' ###########################################################################'
|
2015-03-25 21:36:19 +05:30
|
|
|
close(statUnit)
|
|
|
|
|
2019-05-15 02:14:38 +05:30
|
|
|
call quit(0) ! no complains ;)
|
2012-07-20 21:03:13 +05:30
|
|
|
|
2016-01-27 22:18:27 +05:30
|
|
|
end program DAMASK_spectral
|