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-02-02 17:06:09 +05:30
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
2017-10-05 20:05:34 +05:30
use , intrinsic :: iso_fortran_env , only : &
compiler_version , &
compiler_options
#endif
2018-05-17 15:34:21 +05:30
#include <petsc/finclude/petscsys.h>
use PETScsys
2013-02-23 05:54:30 +05:30
use prec , only : &
pInt , &
2015-08-21 23:21:05 +05:30
pLongInt , &
2013-09-14 16:29:35 +05:30
pReal , &
2016-03-27 00:25:44 +05:30
tol_math_check , &
dNeq
2012-07-20 21:03:13 +05:30
use DAMASK_interface , only : &
DAMASK_interface_init , &
loadCaseFile , &
geometryFile , &
getSolverJobName , &
2018-08-21 02:44:34 +05:30
interface_restartInc
2012-07-19 22:54:56 +05:30
use IO , only : &
2012-07-20 21:03:13 +05:30
IO_isBlank , &
IO_stringPos , &
IO_stringValue , &
IO_floatValue , &
IO_intValue , &
IO_error , &
IO_lc , &
2013-01-10 19:03:43 +05:30
IO_intOut , &
2013-05-08 21:22:29 +05:30
IO_warning , &
2018-08-31 12:28:13 +05:30
IO_timeStamp
2013-05-08 21:22:29 +05:30
use debug , only : &
debug_level , &
debug_spectral , &
debug_levelBasic
2012-10-24 17:01:40 +05:30
use math ! need to include the whole module for FFTW
2015-03-25 21:36:19 +05:30
use mesh , only : &
2015-09-11 14:22:03 +05:30
grid , &
geomSize
2016-01-17 20:33:54 +05:30
use CPFEM2 , only : &
2012-07-20 21:03:13 +05:30
CPFEM_initAll
2012-07-19 22:54:56 +05:30
use FEsolving , only : &
restartWrite , &
restartInc
2012-07-20 21:03:13 +05:30
use numerics , only : &
2015-03-25 21:36:19 +05:30
worldrank , &
worldsize , &
2015-06-03 23:00:31 +05:30
stagItMax , &
2012-10-02 20:56:56 +05:30
maxCutBack , &
2014-09-04 01:29:47 +05:30
spectral_solver , &
2014-03-31 15:34:11 +05:30
continueCalculation
2012-07-19 22:54:56 +05:30
use homogenization , only : &
materialpoint_sizeResults , &
2016-01-17 20:33:54 +05:30
materialpoint_results , &
materialpoint_postResults
2015-06-03 23:00:31 +05:30
use material , only : &
thermal_type , &
damage_type , &
THERMAL_conduction_ID , &
DAMAGE_nonlocal_ID
2015-12-16 02:15:54 +05:30
use spectral_utilities , only : &
2015-06-03 23:00:31 +05:30
utilities_init , &
2012-10-02 20:56:56 +05:30
tSolutionState , &
2015-06-03 23:00:31 +05:30
tLoadCase , &
cutBack , &
nActiveFields , &
FIELD_UNDEFINED_ID , &
FIELD_MECH_ID , &
FIELD_THERMAL_ID , &
2018-02-16 20:06:18 +05:30
FIELD_DAMAGE_ID
2016-01-25 02:35:36 +05:30
use spectral_mech_Basic
use spectral_mech_Polarisation
2015-06-03 23:00:31 +05:30
use spectral_damage
use spectral_thermal
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
2015-08-28 13:08:48 +05:30
integer ( pInt ) , allocatable , dimension ( : ) :: chunkPos
2012-07-20 21:03:13 +05:30
integer ( pInt ) :: &
2016-08-22 19:57:49 +05:30
N_t = 0_pInt , & !< # of time indicators found in load case file
2014-09-04 01:29:47 +05:30
N_n = 0_pInt , & !< # of increment specifiers found in load case file
2013-12-18 14:39:32 +05:30
N_def = 0_pInt !< # of rate of deformation specifiers found in load case file
2013-06-27 00:49:00 +05:30
character ( len = 65536 ) :: &
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
2013-01-08 03:12:00 +05:30
integer ( pInt ) , parameter :: &
subStepFactor = 2_pInt !< 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
stagIterate
2012-12-15 21:51:10 +05:30
integer ( pInt ) :: &
2015-06-03 23:00:31 +05:30
i , j , k , l , field , &
2018-08-31 13:09:33 +05:30
errorID = 0_pInt , &
2012-12-15 21:51:10 +05:30
cutBackLevel = 0_pInt , & !< cut back level \f$ t = \frac{t_{inc}}{2^l} \f$
stepFraction = 0_pInt !< fraction of current time interval
integer ( pInt ) :: &
currentLoadcase = 0_pInt , & !< current load case
inc , & !< current increment in current load case
2014-09-04 01:29:47 +05:30
totalIncsCounter = 0_pInt , & !< total # of increments
convergedCounter = 0_pInt , & !< # of converged increments
notConvergedCounter = 0_pInt , & !< # of non-converged increments
2018-08-31 13:44:33 +05:30
fileUnit = 0_pInt , & !< file unit for reading load case and writing results
myStat , &
2012-12-15 21:51:10 +05:30
statUnit = 0_pInt , & !< file unit for statistics output
2015-08-21 23:21:05 +05:30
lastRestartWritten = 0_pInt , & !< total increment # at which last restart information was written
stagIter
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
2015-08-21 23:21:05 +05:30
integer ( MPI_OFFSET_KIND ) :: fileOffset
integer ( MPI_OFFSET_KIND ) , dimension ( : ) , allocatable :: outputSize
integer ( pInt ) , parameter :: maxByteOut = 2147483647 - 4096 !< limit of one file output write https://trac.mpich.org/projects/mpich/ticket/1742
2016-06-27 21:16:34 +05:30
integer ( pInt ) , parameter :: maxRealOut = maxByteOut / pReal
2015-08-21 23:21:05 +05:30
integer ( pLongInt ) , dimension ( 2 ) :: outputIndex
2018-09-28 11:19:52 +05:30
PetscErrorCode :: ierr
2018-08-31 11:50:23 +05:30
procedure ( basic_init ) , pointer :: &
mech_init
procedure ( basic_forward ) , pointer :: &
mech_forward
procedure ( basic_solution ) , pointer :: &
mech_solution
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
2016-07-29 20:04:51 +05:30
write ( 6 , '(/,a)' ) ' <<<+- DAMASK_spectral init -+>>>'
2018-05-28 10:31:50 +05:30
write ( 6 , '(/,a,/)' ) ' Roters et al., Computational Materials Science, 2018'
2016-07-29 20:04:51 +05:30
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-07-19 22:54:56 +05:30
#include "compilation_info.f90"
2016-08-22 19:57:49 +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
select case ( spectral_solver )
case ( DAMASK_spectral_SolverBasic_label )
mech_init = > basic_init
mech_forward = > basic_forward
mech_solution = > basic_solution
case ( DAMASK_spectral_SolverPolarisation_label )
if ( iand ( debug_level ( debug_spectral ) , debug_levelBasic ) / = 0 ) &
call IO_warning ( 42_pInt , ext_msg = 'debug Divergence' )
mech_init = > polarisation_init
mech_forward = > polarisation_forward
mech_solution = > polarisation_solution
case default
call IO_error ( error_ID = 891_pInt , ext_msg = trim ( spectral_solver ) )
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' )
if ( myStat / = 0_pInt ) call IO_error ( 100_pInt , 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
if ( myStat / = 0_pInt ) 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
currentLoadCase = currentLoadCase + 1_pInt
2015-08-28 13:08:48 +05:30
chunkPos = IO_stringPos ( line )
do i = 1_pInt , chunkPos ( 1 ) ! reading compulsory parameters for loadcase
select case ( IO_lc ( IO_stringValue ( line , chunkPos , i ) ) )
2013-05-13 15:14:23 +05:30
case ( 'l' , 'velocitygrad' , 'velgrad' , 'velocitygradient' , 'fdot' , 'dotf' , 'f' )
N_def = N_def + 1_pInt
2013-02-23 05:54:30 +05:30
case ( 't' , 'time' , 'delta' )
N_t = N_t + 1_pInt
case ( 'n' , 'incs' , 'increments' , 'steps' , 'logincs' , 'logincrements' , 'logsteps' )
N_n = N_n + 1_pInt
end select
2018-08-31 13:49:37 +05:30
enddo
if ( ( N_def / = N_n ) . or . ( N_n / = N_t ) . or . N_n < 1_pInt ) & ! sanity check
call IO_error ( error_ID = 837_pInt , 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
2018-08-31 13:44:33 +05:30
readIn : do i = 1_pInt , chunkPos ( 1 )
2015-08-28 13:08:48 +05:30
select case ( IO_lc ( IO_stringValue ( line , chunkPos , i ) ) )
2013-05-13 15:14:23 +05:30
case ( 'fdot' , 'dotf' , 'l' , 'velocitygrad' , 'velgrad' , 'velocitygradient' , '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
2013-02-06 22:15:34 +05:30
do j = 1_pInt , 9_pInt
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
newLoadCase % deformation % values = math_plain9to33 ( temp_valueVector ) ! values in 3x3 notation
2013-01-24 00:03:46 +05:30
case ( 'p' , 'pk1' , 'piolakirchhoff' , 'stress' , 's' )
2012-07-19 22:54:56 +05:30
temp_valueVector = 0.0_pReal
2013-02-06 22:15:34 +05:30
do j = 1_pInt , 9_pInt
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 )
newLoadCase % stress % values = math_plain9to33 ( temp_valueVector )
2012-07-19 22:54:56 +05:30
case ( 't' , 'time' , 'delta' ) ! increment time
2018-08-31 13:09:33 +05:30
newLoadCase % time = IO_floatValue ( line , chunkPos , i + 1_pInt )
2012-07-19 22:54:56 +05:30
case ( 'n' , 'incs' , 'increments' , 'steps' ) ! number of increments
2018-08-31 13:09:33 +05:30
newLoadCase % incs = IO_intValue ( line , chunkPos , i + 1_pInt )
2012-07-19 22:54:56 +05:30
case ( 'logincs' , 'logincrements' , 'logsteps' ) ! number of increments (switch to log time scaling)
2018-08-31 13:09:33 +05:30
newLoadCase % incs = IO_intValue ( line , chunkPos , i + 1_pInt )
newLoadCase % logscale = 1_pInt
2013-02-23 05:54:30 +05:30
case ( 'freq' , 'frequency' , 'outputfreq' ) ! frequency of result writings
2018-08-31 13:09:33 +05:30
newLoadCase % outputfrequency = IO_intValue ( line , chunkPos , i + 1_pInt )
2012-07-19 22:54:56 +05:30
case ( 'r' , 'restart' , 'restartwrite' ) ! frequency of writing restart information
2018-08-31 13:09:33 +05:30
newLoadCase % restartfrequency = &
2016-08-22 19:57:49 +05:30
max ( 0_pInt , IO_intValue ( line , chunkPos , i + 1_pInt ) )
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
2012-11-28 20:34:05 +05:30
l = 1_pInt ! assuming values given in degrees
k = 1_pInt ! assuming keyword indicating degree/radians present
2015-08-28 13:08:48 +05:30
select case ( IO_lc ( IO_stringValue ( line , chunkPos , i + 1_pInt ) ) )
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
2012-11-28 20:34:05 +05:30
l = 0_pInt
2016-08-22 19:57:49 +05:30
case default
k = 0_pInt
2012-07-19 22:54:56 +05:30
end select
2013-02-06 22:15:34 +05:30
do j = 1_pInt , 3_pInt
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
2012-11-28 20:34:05 +05:30
if ( l == 1_pInt ) temp_valueVector ( 1 : 3 ) = temp_valueVector ( 1 : 3 ) * inRad ! convert to rad
2018-08-31 13:09:33 +05:30
newLoadCase % rotation = math_EulerToR ( temp_valueVector ( 1 : 3 ) ) ! convert rad Eulers to rotation matrix
case ( 'rotation' , 'rot' ) ! assign values for the rotation matrix
2012-07-19 22:54:56 +05:30
temp_valueVector = 0.0_pReal
2013-02-06 22:15:34 +05:30
do j = 1_pInt , 9_pInt
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
2018-08-31 13:09:33 +05:30
newLoadCase % rotation = math_plain9to33 ( temp_valueVector )
2012-07-19 22:54:56 +05:30
end select
2018-08-31 13:44:33 +05:30
enddo readIn
newLoadCase % followFormerTrajectory = merge ( . true . , . false . , currentLoadCase > 1_pInt ) ! 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
write ( 6 , '(1x,a,i6)' ) '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
2015-03-25 21:36:19 +05:30
do j = 1_pInt , 3_pInt
2018-08-31 13:09:33 +05:30
if ( any ( newLoadCase % deformation % maskLogical ( j , 1 : 3 ) . eqv . . true . ) . and . &
any ( newLoadCase % deformation % maskLogical ( j , 1 : 3 ) . eqv . . false . ) ) errorID = 832_pInt ! 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
2015-10-26 22:41:36 +05:30
do i = 1_pInt , 3_pInt ; do j = 1_pInt , 3_pInt
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 . &
newLoadCase % deformation % maskLogical ) ) errorID = 831_pInt ! exclusive or masking only
if ( any ( newLoadCase % stress % maskLogical . and . &
transpose ( newLoadCase % stress % maskLogical ) . and . &
2015-03-25 21:36:19 +05:30
reshape ( [ . false . , . true . , . true . , . true . , . false . , . true . , . true . , . true . , . false . ] , [ 3 , 3 ] ) ) ) &
errorID = 838_pInt ! no rotation is allowed by stress BC
2015-10-26 22:41:36 +05:30
write ( 6 , '(2x,a)' ) 'stress / GPa:'
do i = 1_pInt , 3_pInt ; do j = 1_pInt , 3_pInt
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
2018-08-31 13:09:33 +05:30
if ( any ( abs ( math_mul33x33 ( newLoadCase % rotation , &
transpose ( newLoadCase % rotation ) ) - math_I3 ) > &
2015-03-25 21:36:19 +05:30
reshape ( spread ( tol_math_check , 1 , 9 ) , [ 3 , 3 ] ) ) &
2018-08-31 13:09:33 +05:30
. or . abs ( math_det33 ( newLoadCase % rotation ) ) > &
2015-03-25 21:36:19 +05:30
1.0_pReal + tol_math_check ) errorID = 846_pInt ! given rotation matrix contains strain
2018-08-31 13:09:33 +05:30
if ( any ( dNeq ( newLoadCase % rotation , 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:' , &
2018-08-31 13:09:33 +05:30
transpose ( newLoadCase % rotation )
if ( newLoadCase % time < 0.0_pReal ) errorID = 834_pInt ! negative time increment
write ( 6 , '(2x,a,f12.6)' ) 'time: ' , newLoadCase % time
if ( newLoadCase % incs < 1_pInt ) errorID = 835_pInt ! non-positive incs count
write ( 6 , '(2x,a,i5)' ) 'increments: ' , newLoadCase % incs
if ( newLoadCase % outputfrequency < 1_pInt ) errorID = 836_pInt ! non-positive result frequency
write ( 6 , '(2x,a,i5)' ) 'output frequency: ' , newLoadCase % outputfrequency
write ( 6 , '(2x,a,i5,/)' ) 'restart frequency: ' , newLoadCase % restartfrequency
2015-03-25 21:36:19 +05:30
if ( errorID > 0_pInt ) 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
2015-06-03 23:00:31 +05:30
call Utilities_init ( )
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 )
2015-07-24 20:27:29 +05:30
call spectral_thermal_init
2016-08-22 19:57:49 +05:30
2015-06-03 23:00:31 +05:30
case ( FIELD_DAMAGE_ID )
2018-08-31 11:50:23 +05:30
call spectral_damage_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
2018-08-21 02:44:34 +05:30
writeHeader : if ( interface_restartInc < 1_pInt ) then
2018-08-31 13:44:33 +05:30
open ( newunit = fileUnit , file = trim ( getSolverJobName ( ) ) / / &
2015-03-25 21:36:19 +05:30
'.spectralOut' , form = 'UNFORMATTED' , status = 'REPLACE' )
2018-08-31 13:44:33 +05:30
write ( fileUnit ) 'load:' , trim ( loadCaseFile ) ! ... and write header
2018-09-28 11:19:52 +05:30
write ( fileUnit ) 'workingdir:' , 'n/a'
2018-08-31 13:44:33 +05:30
write ( fileUnit ) 'geometry:' , trim ( geometryFile )
write ( fileUnit ) 'grid:' , grid
write ( fileUnit ) 'size:' , geomSize
write ( fileUnit ) 'materialpoint_sizeResults:' , materialpoint_sizeResults
write ( fileUnit ) 'loadcases:' , size ( loadCases )
write ( fileUnit ) 'frequencies:' , loadCases % outputfrequency ! one entry per LoadCase
write ( fileUnit ) 'times:' , loadCases % time ! one entry per LoadCase
write ( fileUnit ) 'logscales:' , loadCases % logscale
write ( fileUnit ) 'increments:' , loadCases % incs ! one entry per LoadCase
write ( fileUnit ) 'startingIncrement:' , restartInc ! start with writing out the previous inc
write ( fileUnit ) 'eoh'
close ( fileUnit ) ! end of header
2018-07-10 11:54:45 +05:30
open ( newunit = statUnit , file = trim ( getSolverJobName ( ) ) / / &
2015-03-25 21:36:19 +05:30
'.sta' , form = 'FORMATTED' , status = 'REPLACE' )
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 ) &
write ( 6 , '(/,a)' ) ' header of result and statistics file written out'
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
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2015-08-21 23:21:05 +05:30
! prepare MPI parallel out (including opening of file)
allocate ( outputSize ( worldsize ) , source = 0_MPI_OFFSET_KIND )
2016-03-27 12:45:47 +05:30
outputSize ( worldrank + 1 ) = size ( materialpoint_results , kind = MPI_OFFSET_KIND ) * int ( pReal , MPI_OFFSET_KIND )
2016-04-12 14:35:01 +05:30
call MPI_allreduce ( MPI_IN_PLACE , outputSize , worldsize , MPI_LONG , MPI_SUM , PETSC_COMM_WORLD , ierr ) ! get total output size over each process
2016-08-25 21:27:19 +05:30
if ( ierr / = 0_pInt ) call IO_error ( error_ID = 894_pInt , ext_msg = 'MPI_allreduce' )
2018-07-10 11:54:45 +05:30
call MPI_file_open ( PETSC_COMM_WORLD , trim ( getSolverJobName ( ) ) / / '.spectralOut' , &
2015-08-21 23:21:05 +05:30
MPI_MODE_WRONLY + MPI_MODE_APPEND , &
MPI_INFO_NULL , &
2018-08-31 13:44:33 +05:30
fileUnit , &
2015-08-21 23:21:05 +05:30
ierr )
2016-08-25 21:27:19 +05:30
if ( ierr / = 0_pInt ) call IO_error ( error_ID = 894_pInt , ext_msg = 'MPI_file_open' )
2018-08-31 13:44:33 +05:30
call MPI_file_get_position ( fileUnit , fileOffset , ierr ) ! get offset from header
2016-08-25 21:27:19 +05:30
if ( ierr / = 0_pInt ) call IO_error ( error_ID = 894_pInt , ext_msg = 'MPI_file_get_position' )
2015-08-21 23:21:05 +05:30
fileOffset = fileOffset + sum ( outputSize ( 1 : worldrank ) ) ! offset of my process in file (header + processes before me)
2018-08-31 13:44:33 +05:30
call MPI_file_seek ( fileUnit , fileOffset , MPI_SEEK_SET , ierr )
2016-08-25 21:27:19 +05:30
if ( ierr / = 0_pInt ) call IO_error ( error_ID = 894_pInt , ext_msg = 'MPI_file_seek' )
2015-08-21 23:21:05 +05:30
2018-08-21 02:44:34 +05:30
writeUndeformed : if ( interface_restartInc < 1_pInt ) then
2018-02-16 20:06:18 +05:30
write ( 6 , '(1/,a)' ) ' ... writing initial configuration to file ........................'
2016-08-25 21:27:19 +05:30
do i = 1 , size ( materialpoint_results , 3 ) / ( maxByteOut / ( materialpoint_sizeResults * pReal ) ) + 1 ! slice the output of my process in chunks not exceeding the limit for one output
2018-02-16 20:06:18 +05:30
outputIndex = int ( [ ( i - 1_pInt ) * ( ( maxRealOut ) / materialpoint_sizeResults ) + 1_pInt , & ! QUESTION: why not starting i at 0 instead of murky 1?
2016-08-25 21:27:19 +05:30
min ( i * ( ( maxRealOut ) / materialpoint_sizeResults ) , size ( materialpoint_results , 3 ) ) ] , pLongInt )
2018-08-31 13:44:33 +05:30
call MPI_file_write ( fileUnit , reshape ( materialpoint_results ( : , : , outputIndex ( 1 ) : outputIndex ( 2 ) ) , &
2018-05-17 15:34:21 +05:30
[ ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * int ( materialpoint_sizeResults , pLongInt ) ] ) , &
int ( ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * int ( materialpoint_sizeResults , pLongInt ) ) , &
2015-08-21 23:21:05 +05:30
MPI_DOUBLE , MPI_STATUS_IGNORE , ierr )
2016-08-25 21:27:19 +05:30
if ( ierr / = 0_pInt ) call IO_error ( error_ID = 894_pInt , ext_msg = 'MPI_file_write' )
2015-08-21 23:21:05 +05:30
enddo
2016-03-15 03:00:55 +05:30
fileOffset = fileOffset + sum ( outputSize ) ! forward to current file position
2018-08-21 02:44:34 +05:30
endif writeUndeformed
2018-09-28 11:19:52 +05:30
2012-08-03 14:55:48 +05:30
loadCaseLooping : do currentLoadCase = 1_pInt , 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
2012-08-03 14:55:48 +05:30
incLooping : do inc = 1_pInt , loadCases ( currentLoadCase ) % incs
2014-09-04 01:29:47 +05:30
totalIncsCounter = totalIncsCounter + 1_pInt
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
2012-10-02 20:56:56 +05:30
if ( loadCases ( currentLoadCase ) % logscale == 0_pInt ) 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
2018-08-31 13:09:33 +05:30
if ( currentLoadCase == 1_pInt ) then ! 1st load case of logarithmic scale
if ( inc == 1_pInt ) then ! 1st inc of 1st load case of logarithmic scale
2016-08-22 19:57:49 +05:30
timeinc = loadCases ( 1 ) % time * ( 2.0_pReal ** real ( 1_pInt - 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
2012-08-03 14:55:48 +05:30
timeinc = loadCases ( 1 ) % time * ( 2.0_pReal ** real ( inc - 1_pInt - 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 ) ) &
2018-02-17 03:11:07 +05:30
- ( 1.0_pReal + loadCases ( currentLoadCase ) % time / time0 ) ** ( real ( inc - 1_pInt , 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
2018-02-17 03:11:07 +05:30
skipping : if ( totalIncsCounter < = 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
stepFraction = 0_pInt ! 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
stepFraction = stepFraction + 1_pInt ! 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)' ) ' ###########################################################################'
write ( 6 , '(1x,a,es12.5' / / &
2018-02-16 20:06:18 +05:30
',a,' / / IO_intOut ( inc ) / / ',a,' / / IO_intOut ( loadCases ( currentLoadCase ) % incs ) / / &
',a,' / / IO_intOut ( stepFraction ) / / ',a,' / / IO_intOut ( subStepFactor ** cutBackLevel ) / / &
2017-05-24 21:42:36 +05:30
',a,' / / IO_intOut ( currentLoadCase ) / / ',a,' / / IO_intOut ( size ( loadCases ) ) / / ')' ) &
'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 )
2018-02-16 20:06:18 +05:30
write ( incInfo , &
'(a,' / / IO_intOut ( totalIncsCounter ) / / &
',a,' / / IO_intOut ( sum ( loadCases % incs ) ) / / &
',a,' / / IO_intOut ( stepFraction ) / / &
',a,' / / IO_intOut ( subStepFactor ** cutBackLevel ) / / ')' ) &
'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 ( &
2015-06-03 23:00:31 +05:30
guess , timeinc , timeIncOld , remainingLoadCaseTime , &
2016-10-26 02:24:32 +05:30
deformation_BC = loadCases ( currentLoadCase ) % deformation , &
stress_BC = loadCases ( currentLoadCase ) % stress , &
2015-06-03 23:00:31 +05:30
rotation_BC = loadCases ( currentLoadCase ) % rotation )
2016-08-22 19:57:49 +05:30
2018-08-31 11:50:23 +05:30
case ( FIELD_THERMAL_ID ) ; call spectral_thermal_forward ( )
case ( FIELD_DAMAGE_ID ) ; call spectral_damage_forward ( )
2015-06-03 23:00:31 +05:30
end select
2016-08-22 19:57:49 +05:30
enddo
2012-12-15 21:51:10 +05:30
!--------------------------------------------------------------------------------------------------
2015-06-03 23:00:31 +05:30
! solve fields
stagIter = 0_pInt
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 , &
stress_BC = loadCases ( currentLoadCase ) % stress , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
2016-08-22 19:57:49 +05:30
2015-06-03 23:00:31 +05:30
case ( FIELD_THERMAL_ID )
2016-10-22 17:10:45 +05:30
solres ( field ) = spectral_thermal_solution ( timeinc , timeIncOld , remainingLoadCaseTime )
2016-08-22 19:57:49 +05:30
2015-06-03 23:00:31 +05:30
case ( FIELD_DAMAGE_ID )
2016-10-22 17:10:45 +05:30
solres ( field ) = spectral_damage_solution ( timeinc , timeIncOld , remainingLoadCaseTime )
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
stagIter = stagIter + 1_pInt
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
if ( ( continueCalculation . or . all ( solres ( : ) % converged . and . solres ( : ) % stagConverged ) ) & ! don't care or did converge
. and . . not . solres ( 1 ) % termIll ) then ! and acceptable solution found
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 .
stepFraction = ( stepFraction - 1_pInt ) * subStepFactor ! adjust to new denominator
cutBackLevel = cutBackLevel + 1_pInt
time = time - timeinc ! rewind time
timeinc = timeinc / real ( subStepFactor , pReal ) ! cut timestep
write ( 6 , '(/,a)' ) ' cutting back '
else ! no more options to continue
call IO_warning ( 850_pInt )
2018-08-31 13:44:33 +05:30
call MPI_file_close ( fileUnit , ierr )
2018-02-16 20:06:18 +05:30
close ( statUnit )
call quit ( - 1_pInt * ( lastRestartWritten + 1_pInt ) ) ! quit and provide information about last restart inc written
2016-08-22 19:57:49 +05:30
endif
2018-02-16 20:06:18 +05:30
enddo subStepLooping
2012-12-15 21:51:10 +05:30
cutBackLevel = max ( 0_pInt , cutBackLevel - 1_pInt ) ! try half number of subincs next inc
2018-02-16 20:06:18 +05:30
if ( all ( solres ( : ) % converged ) ) then
2012-07-19 22:54:56 +05:30
convergedCounter = convergedCounter + 1_pInt
2018-02-16 20:06:18 +05:30
write ( 6 , '(/,a,' / / IO_intOut ( totalIncsCounter ) / / ',a)' ) & ! report converged inc
2017-05-24 21:42:36 +05:30
' increment ' , totalIncsCounter , ' converged'
2012-07-20 21:03:13 +05:30
else
notConvergedCounter = notConvergedCounter + 1_pInt
2018-02-16 20:06:18 +05:30
write ( 6 , '(/,a,' / / IO_intOut ( totalIncsCounter ) / / ',a)' ) & ! report non-converged inc
' increment ' , totalIncsCounter , ' NOT converged'
2013-02-23 05:54:30 +05:30
endif ; flush ( 6 )
2018-02-16 20:06:18 +05:30
2012-10-02 20:56:56 +05:30
if ( mod ( inc , loadCases ( currentLoadCase ) % outputFrequency ) == 0_pInt ) then ! at output frequency
2018-02-16 20:06:18 +05:30
write ( 6 , '(1/,a)' ) ' ... writing results to file ......................................'
flush ( 6 )
2016-01-17 20:33:54 +05:30
call materialpoint_postResults ( )
2018-08-31 13:44:33 +05:30
call MPI_file_seek ( fileUnit , fileOffset , MPI_SEEK_SET , ierr )
2018-02-16 20:06:18 +05:30
if ( ierr / = 0_pInt ) call IO_error ( 894_pInt , ext_msg = 'MPI_file_seek' )
2015-08-21 23:21:05 +05:30
do i = 1 , size ( materialpoint_results , 3 ) / ( maxByteOut / ( materialpoint_sizeResults * pReal ) ) + 1 ! slice the output of my process in chunks not exceeding the limit for one output
2016-06-27 21:16:34 +05:30
outputIndex = int ( [ ( i - 1_pInt ) * ( ( maxRealOut ) / materialpoint_sizeResults ) + 1_pInt , &
min ( i * ( ( maxRealOut ) / materialpoint_sizeResults ) , size ( materialpoint_results , 3 ) ) ] , pLongInt )
2018-08-31 13:44:33 +05:30
call MPI_file_write ( fileUnit , reshape ( materialpoint_results ( : , : , outputIndex ( 1 ) : outputIndex ( 2 ) ) , &
2018-05-17 15:34:21 +05:30
[ ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * int ( materialpoint_sizeResults , pLongInt ) ] ) , &
int ( ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * int ( materialpoint_sizeResults , pLongInt ) ) , &
2015-08-21 23:21:05 +05:30
MPI_DOUBLE , MPI_STATUS_IGNORE , ierr )
2016-03-27 00:25:44 +05:30
if ( ierr / = 0_pInt ) call IO_error ( 894_pInt , ext_msg = 'MPI_file_write' )
2015-08-21 23:21:05 +05:30
enddo
2016-03-15 03:00:55 +05:30
fileOffset = fileOffset + sum ( outputSize ) ! forward to current file position
2012-07-19 22:54:56 +05:30
endif
2018-02-16 20:06:18 +05:30
if ( loadCases ( currentLoadCase ) % restartFrequency > 0_pInt & ! writing of restart info requested ...
. and . mod ( inc , loadCases ( currentLoadCase ) % restartFrequency ) == 0_pInt ) then ! ... and at frequency of writing restart information
restartWrite = . true . ! set restart parameter for FEsolving
lastRestartWritten = inc ! QUESTION: first call to CPFEM_general will write?
2016-08-22 19:57:49 +05:30
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)' ) ' ###########################################################################'
2018-02-16 20:06:18 +05:30
write ( 6 , '(1x,' / / IO_intOut ( convergedCounter ) / / ',a,' / / IO_intOut ( notConvergedCounter + convergedCounter ) / / ',a,f5.1,a)' ) &
convergedCounter , ' out of ' , &
notConvergedCounter + convergedCounter , ' (' , &
real ( convergedCounter , pReal ) / &
2018-08-31 13:09:33 +05:30
real ( notConvergedCounter + convergedCounter , pReal ) * 10 0.0_pReal , ' %) increments converged!'
2018-02-16 20:06:18 +05:30
flush ( 6 )
2018-08-31 13:44:33 +05:30
call MPI_file_close ( fileUnit , ierr )
2015-03-25 21:36:19 +05:30
close ( statUnit )
2018-09-27 11:53:30 +05:30
if ( notConvergedCounter > 0_pInt ) call quit ( 2_pInt ) ! error if some are not converged
2012-12-15 21:51:10 +05:30
call quit ( 0_pInt ) ! no complains ;)
2012-07-20 21:03:13 +05:30
2016-01-27 22:18:27 +05:30
end program DAMASK_spectral