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
2012-12-15 21:51:10 +05:30
use , intrinsic :: &
iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
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 , &
tol_math_check
2012-07-20 21:03:13 +05:30
use DAMASK_interface , only : &
DAMASK_interface_init , &
loadCaseFile , &
geometryFile , &
getSolverWorkingDirectoryName , &
getSolverJobName , &
appendToOutFile
2012-07-19 22:54:56 +05:30
use IO , only : &
2013-06-27 00:49:00 +05:30
IO_read , &
2012-07-20 21:03:13 +05:30
IO_isBlank , &
IO_open_file , &
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 , &
2013-12-16 17:28:03 +05:30
IO_timeStamp , &
IO_EOF
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 , &
utilities_destroy , &
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 , &
FIELD_DAMAGE_ID
2016-01-25 02:35:36 +05:30
use spectral_mech_Basic
use spectral_mech_AL
use spectral_mech_Polarisation
2015-06-03 23:00:31 +05:30
use spectral_damage
use spectral_thermal
2013-02-20 20:07:12 +05:30
2012-07-19 22:54:56 +05:30
implicit none
2015-03-25 21:36:19 +05:30
2015-08-04 20:34:53 +05:30
#include <petsc/finclude/petscsys.h>
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 ) , parameter :: FILEUNIT = 234_pInt !< file unit, DAMASK IO does not support newunit feature
integer ( pInt ) , allocatable , dimension ( : ) :: chunkPos
2012-07-20 21:03:13 +05:30
integer ( pInt ) :: &
2014-09-04 01:29:47 +05:30
N_t = 0_pInt , & !< # of time indicators found in load case file
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 , &
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 , &
2012-12-15 21:51:10 +05:30
errorID , &
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
2012-12-15 21:51:10 +05:30
resUnit = 0_pInt , & !< file unit for results writing
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
2012-12-15 21:51:10 +05:30
character ( len = 1024 ) :: incInfo !< string parsed to solution with information about current load case
type ( tLoadCase ) , allocatable , dimension ( : ) :: loadCases !< array of all load cases
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
integer ( pLongInt ) , dimension ( 2 ) :: outputIndex
2015-03-25 21:36:19 +05:30
PetscErrorCode :: ierr
2015-08-21 23:21:05 +05:30
external :: &
quit , &
MPI_file_open , &
MPI_file_close , &
MPI_file_seek , &
MPI_file_get_position , &
MPI_file_write , &
MPI_allreduce
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
! init DAMASK (all modules)
2015-07-24 20:27:29 +05:30
call CPFEM_initAll ( el = 1_pInt , ip = 1_pInt )
2015-03-25 21:36:19 +05:30
mainProcess : if ( worldrank == 0 ) then
2016-01-27 22:18:27 +05:30
write ( 6 , '(/,a)' ) ' <<<+- DAMASK_spectral init -+>>>'
2015-03-25 21:36:19 +05:30
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-07-19 22:54:56 +05:30
#include "compilation_info.f90"
2015-03-25 21:36:19 +05:30
endif mainProcess
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 ) )
2012-10-24 17:01:40 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
! reading basic information from load case file and allocate data structure containing load cases
2013-12-12 22:39:59 +05:30
call IO_open_file ( FILEUNIT , trim ( loadCaseFile ) )
rewind ( FILEUNIT )
2012-07-19 22:54:56 +05:30
do
2013-12-12 22:39:59 +05:30
line = IO_read ( FILEUNIT )
2013-12-16 17:28:03 +05:30
if ( trim ( line ) == IO_EOF ) exit
2012-07-19 22:54:56 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
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
2015-10-13 00:32:42 +05:30
enddo ! count all identifiers to allocate memory and do sanity check
2012-07-19 22:54:56 +05:30
enddo
2015-10-13 00:32:42 +05:30
if ( ( N_def / = N_n ) . or . ( N_n / = N_t ) . or . N_n < 1_pInt ) & ! sanity check
2013-06-27 00:49:00 +05:30
call IO_error ( error_ID = 837_pInt , ext_msg = trim ( loadCaseFile ) ) ! error message for incomplete loadcase
2012-12-15 21:51:10 +05:30
allocate ( loadCases ( N_n ) ) ! array of load cases
2012-08-03 14:55:48 +05:30
loadCases % P % myType = 'p'
2015-06-03 23:00:31 +05:30
do i = 1 , size ( loadCases )
allocate ( loadCases ( i ) % ID ( nActiveFields ) )
field = 1
loadCases ( i ) % ID ( field ) = FIELD_MECH_ID ! mechanical active by default
if ( any ( thermal_type == THERMAL_conduction_ID ) ) then ! thermal field active
field = field + 1
loadCases ( i ) % ID ( field ) = FIELD_THERMAL_ID
endif
if ( any ( damage_type == DAMAGE_nonlocal_ID ) ) then ! damage field active
field = field + 1
loadCases ( i ) % ID ( field ) = FIELD_DAMAGE_ID
endif
enddo
2012-07-19 22:54:56 +05:30
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-19 22:54:56 +05:30
! reading the load case and assign values to the allocated data structure
2013-12-12 22:39:59 +05:30
rewind ( FILEUNIT )
2012-07-19 22:54:56 +05:30
do
2013-12-12 22:39:59 +05:30
line = IO_read ( FILEUNIT )
2013-12-16 17:28:03 +05:30
if ( trim ( line ) == IO_EOF ) exit
2012-07-19 22:54:56 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2012-08-03 14:55:48 +05:30
currentLoadCase = currentLoadCase + 1_pInt
2015-08-28 13:08:48 +05:30
chunkPos = IO_stringPos ( line )
do i = 1_pInt , chunkPos ( 1 )
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
2012-08-03 14:55:48 +05:30
loadCases ( currentLoadCase ) % deformation % myType = 'fdot'
2015-08-28 13:08:48 +05:30
else if ( IO_lc ( IO_stringValue ( line , chunkPos , i ) ) == 'f' ) then
2013-05-13 15:14:23 +05:30
loadCases ( currentLoadCase ) % deformation % myType = 'f'
2012-11-29 18:56:17 +05:30
else
loadCases ( currentLoadCase ) % 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 *
2013-02-06 22:15:34 +05:30
enddo
2012-12-15 21:51:10 +05:30
do j = 1_pInt , 9_pInt
2015-08-28 13:08:48 +05:30
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
2012-12-15 21:51:10 +05:30
loadCases ( currentLoadCase ) % deformation % maskLogical = & ! logical mask in 3x3 notation
transpose ( reshape ( temp_maskVector , [ 3 , 3 ] ) )
loadCases ( currentLoadCase ) % deformation % maskFloat = & ! float (1.0/0.0) mask in 3x3 notation
merge ( ones , zeros , loadCases ( currentLoadCase ) % deformation % maskLogical )
loadCases ( currentLoadCase ) % 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
2013-02-06 22:15:34 +05:30
enddo
2012-07-20 21:03:13 +05:30
do j = 1_pInt , 9_pInt
2015-08-28 13:08:48 +05:30
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
2012-08-03 14:55:48 +05:30
loadCases ( currentLoadCase ) % P % maskLogical = transpose ( reshape ( temp_maskVector , [ 3 , 3 ] ) )
2012-10-19 14:14:21 +05:30
loadCases ( currentLoadCase ) % P % maskFloat = merge ( ones , zeros , &
2012-08-03 14:55:48 +05:30
loadCases ( currentLoadCase ) % P % maskLogical )
loadCases ( currentLoadCase ) % P % values = math_plain9to33 ( temp_valueVector )
2012-07-19 22:54:56 +05:30
case ( 't' , 'time' , 'delta' ) ! increment time
2015-08-28 13:08:48 +05:30
loadCases ( currentLoadCase ) % time = IO_floatValue ( line , chunkPos , i + 1_pInt )
2012-07-19 22:54:56 +05:30
case ( 'n' , 'incs' , 'increments' , 'steps' ) ! number of increments
2015-08-28 13:08:48 +05:30
loadCases ( currentLoadCase ) % 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)
2015-08-28 13:08:48 +05:30
loadCases ( currentLoadCase ) % incs = IO_intValue ( line , chunkPos , i + 1_pInt )
2012-08-03 14:55:48 +05:30
loadCases ( currentLoadCase ) % logscale = 1_pInt
2013-02-23 05:54:30 +05:30
case ( 'freq' , 'frequency' , 'outputfreq' ) ! frequency of result writings
2015-08-28 13:08:48 +05:30
loadCases ( currentLoadCase ) % 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
2012-12-15 21:51:10 +05:30
loadCases ( currentLoadCase ) % restartfrequency = &
2015-08-28 13:08:48 +05:30
max ( 0_pInt , IO_intValue ( line , chunkPos , i + 1_pInt ) )
2012-07-19 22:54:56 +05:30
case ( 'guessreset' , 'dropguessing' )
2012-10-24 17:01:40 +05:30
loadCases ( currentLoadCase ) % followFormerTrajectory = . false . ! do not continue to predict deformation along former trajectory
2012-08-03 14:55:48 +05:30
case ( 'euler' ) ! rotation of currentLoadCase 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' )
2012-11-28 20:34:05 +05:30
case ( 'rad' , 'radian' ) ! don't convert from degree to radian
l = 0_pInt
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
2012-09-06 19:35:28 +05:30
loadCases ( currentLoadCase ) % rotation = math_EulerToR ( temp_valueVector ( 1 : 3 ) ) ! convert rad Eulers to rotation matrix
2012-08-03 14:55:48 +05:30
case ( 'rotation' , 'rot' ) ! assign values for the rotation of currentLoadCase 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
2012-08-03 14:55:48 +05:30
loadCases ( currentLoadCase ) % rotation = math_plain9to33 ( temp_valueVector )
2012-07-19 22:54:56 +05:30
end select
enddo ; enddo
2013-12-18 14:39:32 +05:30
close ( FILEUNIT )
2012-08-03 14:55:48 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! consistency checks and output of load case
2012-10-24 17:01:40 +05:30
loadCases ( 1 ) % followFormerTrajectory = . false . ! cannot guess along trajectory for first inc of first currentLoadCase
2012-07-19 22:54:56 +05:30
errorID = 0_pInt
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) then
checkLoadcases : do currentLoadCase = 1_pInt , size ( loadCases )
write ( loadcase_string , '(i6)' ) currentLoadCase
write ( 6 , '(1x,a,i6)' ) 'load case: ' , currentLoadCase
if ( . not . loadCases ( currentLoadCase ) % followFormerTrajectory ) &
write ( 6 , '(2x,a)' ) 'drop guessing along trajectory'
if ( loadCases ( currentLoadCase ) % deformation % myType == 'l' ) then
do j = 1_pInt , 3_pInt
if ( any ( loadCases ( currentLoadCase ) % deformation % maskLogical ( j , 1 : 3 ) . eqv . . true . ) . and . &
any ( loadCases ( currentLoadCase ) % deformation % maskLogical ( j , 1 : 3 ) . eqv . . false . ) ) &
errorID = 832_pInt ! each row should be either fully or not at all defined
enddo
write ( 6 , '(2x,a)' ) 'velocity gradient:'
else if ( loadCases ( currentLoadCase ) % deformation % myType == 'f' ) then
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
if ( loadCases ( currentLoadCase ) % deformation % maskLogical ( i , j ) ) then
write ( 6 , '(2x,f12.7)' , advance = 'no' ) loadCases ( currentLoadCase ) % deformation % values ( i , j )
else
write ( 6 , '(2x,12a)' , advance = 'no' ) ' * '
endif
enddo ; write ( 6 , '(/)' , advance = 'no' )
enddo
2015-03-25 21:36:19 +05:30
if ( any ( loadCases ( currentLoadCase ) % P % maskLogical . eqv . &
loadCases ( currentLoadCase ) % deformation % maskLogical ) ) errorID = 831_pInt ! exclusive or masking only
if ( any ( loadCases ( currentLoadCase ) % P % maskLogical . and . &
transpose ( loadCases ( currentLoadCase ) % P % maskLogical ) . and . &
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
2016-03-09 22:10:12 +05:30
if ( loadCases ( currentLoadCase ) % P % maskLogical ( i , j ) ) then
2015-10-26 22:41:36 +05:30
write ( 6 , '(2x,f12.7)' , advance = 'no' ) loadCases ( currentLoadCase ) % P % values ( i , j ) * 1e-9_pReal
else
write ( 6 , '(2x,12a)' , advance = 'no' ) ' * '
endif
enddo ; write ( 6 , '(/)' , advance = 'no' )
enddo
2015-03-25 21:36:19 +05:30
if ( any ( abs ( math_mul33x33 ( loadCases ( currentLoadCase ) % rotation , &
math_transpose33 ( loadCases ( currentLoadCase ) % rotation ) ) - math_I3 ) > &
reshape ( spread ( tol_math_check , 1 , 9 ) , [ 3 , 3 ] ) ) &
. or . abs ( math_det33 ( loadCases ( currentLoadCase ) % rotation ) ) > &
1.0_pReal + tol_math_check ) errorID = 846_pInt ! given rotation matrix contains strain
if ( any ( loadCases ( currentLoadCase ) % rotation / = math_I3 ) ) &
write ( 6 , '(2x,a,/,3(3(3x,f12.7,1x)/))' , advance = 'no' ) 'rotation of loadframe:' , &
math_transpose33 ( loadCases ( currentLoadCase ) % rotation )
if ( loadCases ( currentLoadCase ) % time < 0.0_pReal ) errorID = 834_pInt ! negative time increment
write ( 6 , '(2x,a,f12.6)' ) 'time: ' , loadCases ( currentLoadCase ) % time
if ( loadCases ( currentLoadCase ) % incs < 1_pInt ) errorID = 835_pInt ! non-positive incs count
write ( 6 , '(2x,a,i5)' ) 'increments: ' , loadCases ( currentLoadCase ) % incs
if ( loadCases ( currentLoadCase ) % outputfrequency < 1_pInt ) errorID = 836_pInt ! non-positive result frequency
write ( 6 , '(2x,a,i5)' ) 'output frequency: ' , &
loadCases ( currentLoadCase ) % outputfrequency
write ( 6 , '(2x,a,i5,/)' ) 'restart frequency: ' , &
loadCases ( currentLoadCase ) % restartfrequency
if ( errorID > 0_pInt ) call IO_error ( error_ID = errorID , ext_msg = loadcase_string ) ! exit with error message
enddo checkLoadcases
endif
2012-07-20 21:03:13 +05:30
2012-12-15 21:51:10 +05:30
!--------------------------------------------------------------------------------------------------
! doing initialization depending on selected solver
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 )
select case ( spectral_solver )
case ( DAMASK_spectral_SolverBasicPETSc_label )
2015-07-24 20:27:29 +05:30
call basicPETSc_init
2015-06-03 23:00:31 +05:30
case ( DAMASK_spectral_SolverAL_label )
if ( iand ( debug_level ( debug_spectral ) , debug_levelBasic ) / = 0 . and . worldrank == 0_pInt ) &
call IO_warning ( 42_pInt , ext_msg = 'debug Divergence' )
2015-07-24 20:27:29 +05:30
call AL_init
2015-06-03 23:00:31 +05:30
case ( DAMASK_spectral_SolverPolarisation_label )
if ( iand ( debug_level ( debug_spectral ) , debug_levelBasic ) / = 0 . and . worldrank == 0_pInt ) &
call IO_warning ( 42_pInt , ext_msg = 'debug Divergence' )
2015-07-24 20:27:29 +05:30
call Polarisation_init
2015-06-03 23:00:31 +05:30
case default
call IO_error ( error_ID = 891 , ext_msg = trim ( spectral_solver ) )
end select
case ( FIELD_THERMAL_ID )
2015-07-24 20:27:29 +05:30
call spectral_thermal_init
2015-06-03 23:00:31 +05:30
case ( FIELD_DAMAGE_ID )
call spectral_damage_init ( )
end select
enddo
2012-08-29 00:49:47 +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
if ( . not . appendToOutFile ) then ! after restart, append to existing results file
2015-03-25 21:36:19 +05:30
open ( newunit = resUnit , file = trim ( getSolverWorkingDirectoryName ( ) ) / / trim ( getSolverJobName ( ) ) / / &
'.spectralOut' , form = 'UNFORMATTED' , status = 'REPLACE' )
write ( resUnit ) 'load:' , trim ( loadCaseFile ) ! ... and write header
write ( resUnit ) 'workingdir:' , trim ( getSolverWorkingDirectoryName ( ) )
write ( resUnit ) 'geometry:' , trim ( geometryFile )
2015-09-11 14:22:03 +05:30
write ( resUnit ) 'grid:' , grid
write ( resUnit ) 'size:' , geomSize
2015-03-25 21:36:19 +05:30
write ( resUnit ) 'materialpoint_sizeResults:' , materialpoint_sizeResults
write ( resUnit ) 'loadcases:' , size ( loadCases )
2015-04-01 21:34:33 +05:30
write ( resUnit ) 'frequencies:' , loadCases % outputfrequency ! one entry per LoadCase
write ( resUnit ) 'times:' , loadCases % time ! one entry per LoadCase
2015-03-25 21:36:19 +05:30
write ( resUnit ) 'logscales:' , loadCases % logscale
2015-04-01 21:34:33 +05:30
write ( resUnit ) 'increments:' , loadCases % incs ! one entry per LoadCase
2015-03-25 21:36:19 +05:30
write ( resUnit ) 'startingIncrement:' , restartInc - 1_pInt ! start with writing out the previous inc
write ( resUnit ) 'eoh'
close ( resUnit ) ! end of header
open ( newunit = statUnit , file = trim ( getSolverWorkingDirectoryName ( ) ) / / trim ( getSolverJobName ( ) ) / / &
'.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 )
else ! open new files ...
open ( newunit = statUnit , file = trim ( getSolverWorkingDirectoryName ( ) ) / / trim ( getSolverJobName ( ) ) / / &
'.sta' , form = 'FORMATTED' , position = 'APPEND' , status = 'OLD' )
2015-03-25 21:36:19 +05:30
endif
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 )
outputSize ( worldrank + 1 ) = int ( size ( materialpoint_results ) * pReal , MPI_OFFSET_KIND )
call MPI_allreduce ( MPI_IN_PLACE , outputSize , worldsize , MPI_INT , MPI_SUM , PETSC_COMM_WORLD , ierr ) ! get total output size over each process
call MPI_file_open ( PETSC_COMM_WORLD , &
trim ( getSolverWorkingDirectoryName ( ) ) / / trim ( getSolverJobName ( ) ) / / '.spectralOut' , &
MPI_MODE_WRONLY + MPI_MODE_APPEND , &
MPI_INFO_NULL , &
resUnit , &
ierr )
call MPI_file_get_position ( resUnit , fileOffset , ierr ) ! get offset from header
fileOffset = fileOffset + sum ( outputSize ( 1 : worldrank ) ) ! offset of my process in file (header + processes before me)
call MPI_file_seek ( resUnit , fileOffset , MPI_SEEK_SET , ierr )
if ( . not . appendToOutFile ) then ! if not restarting, write 0th increment
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-03-15 03:00:55 +05:30
outputIndex = int ( [ ( i - 1_pInt ) * ( ( maxByteOut / pReal ) / materialpoint_sizeResults ) + 1_pInt , &
min ( i * ( ( maxByteOut / pReal ) / materialpoint_sizeResults ) , size ( materialpoint_results , 3 ) ) ] , pLongInt )
2015-08-21 23:21:05 +05:30
call MPI_file_write ( resUnit , reshape ( materialpoint_results ( : , : , outputIndex ( 1 ) : outputIndex ( 2 ) ) , &
[ ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * materialpoint_sizeResults ] ) , &
( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * materialpoint_sizeResults , &
MPI_DOUBLE , MPI_STATUS_IGNORE , ierr )
enddo
2016-03-15 03:00:55 +05:30
fileOffset = fileOffset + sum ( outputSize ) ! forward to current file position
2015-09-01 22:23:48 +05:30
if ( worldrank == 0 ) &
write ( 6 , '(1/,a)' ) ' ... writing initial configuration to file ........................'
2015-08-21 23:21:05 +05:30
endif
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
! loopping over loadcases
loadCaseLooping : do currentLoadCase = 1_pInt , size ( loadCases )
time0 = time ! currentLoadCase 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
!--------------------------------------------------------------------------------------------------
! loop oper incs defined in input file for current currentLoadCase
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
2013-05-13 15:14:23 +05:30
timeIncOld = timeinc
2012-10-02 20:56:56 +05:30
if ( loadCases ( currentLoadCase ) % logscale == 0_pInt ) then ! linear scale
timeinc = loadCases ( currentLoadCase ) % time / loadCases ( currentLoadCase ) % incs ! only valid for given linear time scale. will be overwritten later in case loglinear scale is used
2012-07-19 22:54:56 +05:30
else
2012-10-02 20:56:56 +05:30
if ( currentLoadCase == 1_pInt ) then ! 1st currentLoadCase of logarithmic scale
2012-08-03 14:55:48 +05:30
if ( inc == 1_pInt ) then ! 1st inc of 1st currentLoadCase of logarithmic scale
2012-10-02 20:56:56 +05:30
timeinc = loadCases ( 1 ) % time * ( 2.0_pReal ** real ( 1_pInt - loadCases ( 1 ) % incs , pReal ) ) ! assume 1st inc is equal to 2nd
2012-08-03 14:55:48 +05:30
else ! not-1st inc of 1st currentLoadCase of logarithmic scale
timeinc = loadCases ( 1 ) % time * ( 2.0_pReal ** real ( inc - 1_pInt - loadCases ( 1 ) % incs , pReal ) )
2012-07-19 22:54:56 +05:30
endif
2012-08-03 14:55:48 +05:30
else ! not-1st currentLoadCase of logarithmic scale
2012-12-15 21:51:10 +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_pInt ) , pReal ) / &
real ( loadCases ( currentLoadCase ) % incs , pReal ) ) )
2012-07-19 22:54:56 +05:30
endif
endif
2012-12-15 21:51:10 +05:30
timeinc = timeinc / 2.0_pReal ** real ( cutBackLevel , pReal ) ! depending on cut back level, decrease time step
2012-07-19 22:54:56 +05:30
2013-01-29 15:58:01 +05:30
forwarding : if ( totalIncsCounter > = restartInc ) then
2012-10-02 20:56:56 +05:30
stepFraction = 0_pInt
2013-02-23 05:54:30 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-12-15 21:51:10 +05:30
! loop over sub incs
2013-01-08 03:12:00 +05:30
subIncLooping : do while ( stepFraction / subStepFactor ** cutBackLevel < 1_pInt )
2012-12-15 21:51:10 +05:30
time = time + timeinc ! forward time
2012-10-02 20:56:56 +05:30
stepFraction = stepFraction + 1_pInt
2013-05-13 15:14:23 +05:30
remainingLoadCaseTime = time0 - time + loadCases ( currentLoadCase ) % time + timeInc
2013-12-18 15:05:05 +05:30
!--------------------------------------------------------------------------------------------------
2012-12-15 21:51:10 +05:30
! report begin of new increment
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) then
write ( 6 , '(/,a)' ) ' ###########################################################################'
write ( 6 , '(1x,a,es12.5' / / &
',a,' / / IO_intOut ( inc ) / / ',a,' / / IO_intOut ( loadCases ( currentLoadCase ) % incs ) / / &
',a,' / / IO_intOut ( stepFraction ) / / ',a,' / / IO_intOut ( subStepFactor ** cutBackLevel ) / / &
',a,' / / IO_intOut ( currentLoadCase ) / / ',a,' / / IO_intOut ( size ( loadCases ) ) / / ')' ) &
'Time' , time , &
's: Increment ' , inc , '/' , loadCases ( currentLoadCase ) % incs , &
'-' , stepFraction , '/' , subStepFactor ** cutBackLevel , &
' of load case ' , currentLoadCase , '/' , size ( loadCases )
flush ( 6 )
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 ) , &
'-' , stepFraction , '/' , subStepFactor ** cutBackLevel
endif
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
! forward fields
do field = 1 , nActiveFields
select case ( loadCases ( currentLoadCase ) % ID ( field ) )
case ( FIELD_MECH_ID )
select case ( spectral_solver )
case ( DAMASK_spectral_SolverBasicPETSc_label )
call BasicPETSc_forward ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
P_BC = loadCases ( currentLoadCase ) % P , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
case ( DAMASK_spectral_SolverAL_label )
call AL_forward ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
P_BC = loadCases ( currentLoadCase ) % P , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
case ( DAMASK_spectral_SolverPolarisation_label )
call Polarisation_forward ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
P_BC = loadCases ( currentLoadCase ) % P , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
end select
case ( FIELD_THERMAL_ID )
call spectral_thermal_forward ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime )
case ( FIELD_DAMAGE_ID )
call spectral_damage_forward ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime )
end select
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 )
select case ( spectral_solver )
case ( DAMASK_spectral_SolverBasicPETSc_label )
solres ( field ) = BasicPETSC_solution ( &
incInfo , guess , timeinc , timeIncOld , remainingLoadCaseTime , &
P_BC = loadCases ( currentLoadCase ) % P , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
case ( DAMASK_spectral_SolverAL_label )
solres ( field ) = AL_solution ( &
incInfo , guess , timeinc , timeIncOld , remainingLoadCaseTime , &
P_BC = loadCases ( currentLoadCase ) % P , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
case ( DAMASK_spectral_SolverPolarisation_label )
solres ( field ) = Polarisation_solution ( &
incInfo , guess , timeinc , timeIncOld , remainingLoadCaseTime , &
P_BC = loadCases ( currentLoadCase ) % P , &
F_BC = loadCases ( currentLoadCase ) % deformation , &
rotation_BC = loadCases ( currentLoadCase ) % rotation )
end select
case ( FIELD_THERMAL_ID )
solres ( field ) = spectral_thermal_solution ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime )
case ( FIELD_DAMAGE_ID )
solres ( field ) = spectral_damage_solution ( &
guess , timeinc , timeIncOld , remainingLoadCaseTime )
end select
if ( . not . solres ( field ) % converged ) exit ! no solution found
enddo
stagIter = stagIter + 1_pInt
stagIterate = stagIter < stagItMax . and . &
all ( solres ( : ) % converged ) . and . &
. not . all ( solres ( : ) % stagConverged )
enddo
2013-02-23 05:54:30 +05:30
2012-12-15 21:51:10 +05:30
!--------------------------------------------------------------------------------------------------
! check solution
2013-04-26 14:43:39 +05:30
cutBack = . False .
2016-01-18 01:12:24 +05:30
if ( solres ( 1 ) % termIll . or . . not . all ( solres ( : ) % converged . and . solres ( : ) % stagConverged ) ) then ! no solution found
2012-10-02 20:56:56 +05:30
if ( cutBackLevel < maxCutBack ) then ! do cut back
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) write ( 6 , '(/,a)' ) ' cut back detected'
2012-10-02 20:56:56 +05:30
cutBack = . True .
2013-01-08 03:12:00 +05:30
stepFraction = ( stepFraction - 1_pInt ) * subStepFactor ! adjust to new denominator
2012-10-02 20:56:56 +05:30
cutBackLevel = cutBackLevel + 1_pInt
2012-12-15 21:51:10 +05:30
time = time - timeinc ! rewind time
2012-10-02 20:56:56 +05:30
timeinc = timeinc / 2.0_pReal
2015-06-03 23:00:31 +05:30
elseif ( solres ( 1 ) % termIll ) then ! material point model cannot find a solution, exit in any casy
2014-03-31 15:34:11 +05:30
call IO_warning ( 850_pInt )
call quit ( - 1_pInt * ( lastRestartWritten + 1_pInt ) ) ! quit and provide information about last restart inc written (e.g. for regridding)
elseif ( continueCalculation == 1_pInt ) then
guess = . true . ! accept non converged BVP solution
else ! default behavior, exit if spectral solver does not converge
call IO_warning ( 850_pInt )
2016-01-18 01:12:24 +05:30
call quit ( - 1_pInt * ( lastRestartWritten + 1_pInt ) ) ! quit and provide information about last restart inc written (e.g. for regridding)
2012-10-02 20:56:56 +05:30
endif
else
2012-11-09 01:02:00 +05:30
guess = . true . ! start guessing after first converged (sub)inc
2012-10-02 20:56:56 +05:30
endif
2015-03-25 21:36:19 +05:30
if ( . not . cutBack ) then
if ( worldrank == 0 ) then
write ( statUnit , * ) totalIncsCounter , time , cutBackLevel , &
solres % converged , solres % iterationsNeeded ! write statistics about accepted solution
flush ( statUnit )
endif
endif
2012-10-02 20:56:56 +05:30
enddo subIncLooping
2012-12-15 21:51:10 +05:30
cutBackLevel = max ( 0_pInt , cutBackLevel - 1_pInt ) ! try half number of subincs next inc
2015-06-03 23:00:31 +05:30
if ( all ( solres ( : ) % converged ) ) then ! report converged inc
2012-07-19 22:54:56 +05:30
convergedCounter = convergedCounter + 1_pInt
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) &
write ( 6 , '(/,a,' / / IO_intOut ( totalIncsCounter ) / / ',a)' ) &
2013-01-10 03:49:32 +05:30
' increment ' , totalIncsCounter , ' converged'
2012-07-20 21:03:13 +05:30
else
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) &
write ( 6 , '(/,a,' / / IO_intOut ( totalIncsCounter ) / / ',a)' ) & ! report non-converged inc
2013-01-10 03:49:32 +05:30
' increment ' , totalIncsCounter , ' NOT converged'
2012-07-20 21:03:13 +05:30
notConvergedCounter = notConvergedCounter + 1_pInt
2013-02-23 05:54:30 +05:30
endif ; flush ( 6 )
2012-10-02 20:56:56 +05:30
if ( mod ( inc , loadCases ( currentLoadCase ) % outputFrequency ) == 0_pInt ) then ! at output frequency
2015-03-25 21:36:19 +05:30
if ( worldrank == 0 ) &
write ( 6 , '(1/,a)' ) ' ... writing results to file ......................................'
2016-01-17 20:33:54 +05:30
call materialpoint_postResults ( )
2015-08-21 23:21:05 +05:30
call MPI_file_seek ( resUnit , fileOffset , MPI_SEEK_SET , ierr )
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-03-15 03:00:55 +05:30
outputIndex = int ( [ ( i - 1_pInt ) * ( ( maxByteOut / pReal ) / materialpoint_sizeResults ) + 1_pInt , &
min ( i * ( ( maxByteOut / pReal ) / materialpoint_sizeResults ) , size ( materialpoint_results , 3 ) ) ] , pLongInt )
2015-08-21 23:21:05 +05:30
call MPI_file_write ( resUnit , reshape ( materialpoint_results ( : , : , outputIndex ( 1 ) : outputIndex ( 2 ) ) , &
[ ( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * materialpoint_sizeResults ] ) , &
( outputIndex ( 2 ) - outputIndex ( 1 ) + 1 ) * materialpoint_sizeResults , &
MPI_DOUBLE , MPI_STATUS_IGNORE , ierr )
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
2013-06-11 22:05:04 +05:30
if ( loadCases ( currentLoadCase ) % restartFrequency > 0_pInt . and . & ! at frequency of writing restart information set restart parameter for FEsolving
2015-03-25 21:36:19 +05:30
mod ( inc , loadCases ( currentLoadCase ) % restartFrequency ) == 0_pInt ) then ! first call to CPFEM_general will write?
2012-10-24 17:01:40 +05:30
restartWrite = . true .
2012-12-14 20:48:04 +05:30
lastRestartWritten = inc
2012-10-24 17:01:40 +05:30
endif
2013-01-29 15:58:01 +05:30
else forwarding
2012-10-24 17:01:40 +05:30
time = time + timeinc
2012-11-09 01:02:00 +05:30
guess = . true .
2013-01-29 15:58:01 +05:30
endif forwarding
2012-07-20 21:03:13 +05:30
2012-07-24 22:37:10 +05:30
enddo incLooping
enddo loadCaseLooping
2015-10-26 22:41:36 +05:30
2015-03-25 21:36:19 +05:30
!--------------------------------------------------------------------------------------------------
! report summary of whole calculation
if ( worldrank == 0 ) then
write ( 6 , '(/,a)' ) ' ###########################################################################'
write ( 6 , '(1x,i6.6,a,i6.6,a,f5.1,a)' ) convergedCounter , ' out of ' , &
notConvergedCounter + convergedCounter , ' (' , &
real ( convergedCounter , pReal ) / &
real ( notConvergedCounter + convergedCounter , pReal ) * 10 0.0_pReal , &
' %) increments converged!'
endif
call MPI_file_close ( resUnit , ierr )
close ( statUnit )
2015-06-03 23:00:31 +05:30
do field = 1 , nActiveFields
select case ( loadCases ( 1 ) % ID ( field ) )
case ( FIELD_MECH_ID )
select case ( spectral_solver )
case ( DAMASK_spectral_SolverBasicPETSc_label )
call BasicPETSC_destroy ( )
case ( DAMASK_spectral_SolverAL_label )
call AL_destroy ( )
case ( DAMASK_spectral_SolverPolarisation_label )
call Polarisation_destroy ( )
end select
case ( FIELD_THERMAL_ID )
call spectral_thermal_destroy ( )
case ( FIELD_DAMAGE_ID )
call spectral_damage_destroy ( )
end select
enddo
call utilities_destroy ( )
2015-10-26 22:41:36 +05:30
2015-03-26 19:13:18 +05:30
call PetscFinalize ( ierr ) ; CHKERRQ ( ierr )
2015-10-26 22:41:36 +05:30
2012-12-15 21:51:10 +05:30
if ( notConvergedCounter > 0_pInt ) call quit ( 3_pInt ) ! error if some are not converged
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
2012-07-20 21:03:13 +05:30
2013-02-23 05:54:30 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief quit subroutine to mimic behavior of FEM solvers
!> @details exits the Spectral solver and reports time and duration. Exit code 0 signals
!> everything went fine. Exit code 1 signals an error, message according to IO_error. Exit code
!> 2 signals request for regridding, increment of last saved restart information is written to
!> stderr. Exit code 3 signals no severe problems, but some increments did not converge
!--------------------------------------------------------------------------------------------------
2013-02-19 20:24:34 +05:30
subroutine quit ( stop_id )
use prec , only : &
pInt
2015-03-25 21:36:19 +05:30
use numerics , only : &
worldrank
2015-10-26 22:41:36 +05:30
2013-02-19 20:24:34 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: stop_id
integer , dimension ( 8 ) :: dateAndTime ! type default integer
2015-03-25 21:36:19 +05:30
if ( worldrank == 0_pInt ) then
call date_and_time ( values = dateAndTime )
write ( 6 , '(/,a)' ) 'DAMASK terminated on:'
write ( 6 , '(a,2(i2.2,a),i4.4)' ) 'Date: ' , dateAndTime ( 3 ) , '/' , &
dateAndTime ( 2 ) , '/' , &
2015-10-26 22:41:36 +05:30
dateAndTime ( 1 )
2015-03-25 21:36:19 +05:30
write ( 6 , '(a,2(i2.2,a),i2.2)' ) 'Time: ' , dateAndTime ( 5 ) , ':' , &
dateAndTime ( 6 ) , ':' , &
2015-10-26 22:41:36 +05:30
dateAndTime ( 7 )
2015-03-25 21:36:19 +05:30
endif
2013-02-19 20:24:34 +05:30
if ( stop_id == 0_pInt ) stop 0 ! normal termination
if ( stop_id < 0_pInt ) then ! trigger regridding
2015-03-25 21:36:19 +05:30
if ( worldrank == 0_pInt ) &
write ( 0 , '(a,i6)' ) 'restart information available at ' , stop_id * ( - 1_pInt )
2013-02-19 20:24:34 +05:30
stop 2
endif
2013-02-23 05:54:30 +05:30
if ( stop_id == 3_pInt ) stop 3 ! not all incs converged
2013-02-19 20:24:34 +05:30
stop 1 ! error (message from IO_error)
2013-02-23 05:54:30 +05:30
end subroutine quit