DAMASK_EICMD/src/grid/DAMASK_grid.f90

595 lines
35 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Driver controlling inner and outer load case looping of the various spectral solvers
!> @details doing cutbacking, forwarding in case of restart, reporting statistics, writing
!> results
!--------------------------------------------------------------------------------------------------
program DAMASK_spectral
#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 mesh
use CPFEM2
use FEsolving
use numerics
use homogenization
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
use HDF5_utilities
2018-11-18 17:11:05 +05:30
use results
2019-09-21 05:22:55 +05:30
use rotations
implicit none
!--------------------------------------------------------------------------------------------------
! variables related to information from load case and geom file
real(pReal), dimension(9) :: temp_valueVector = 0.0_pReal !< temporarily from loadcase file when reading in tensors (initialize to 0.0)
logical, dimension(9) :: temp_maskVector = .false. !< temporarily from loadcase file when reading in tensors
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
character(len=65536) :: &
line
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
real(pReal), dimension(3,3), parameter :: &
ones = 1.0_pReal, &
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
real(pReal) :: &
time = 0.0_pReal, & !< elapsed time
time0 = 0.0_pReal, & !< begin of interval
timeinc = 1.0_pReal, & !< current time interval
timeIncOld = 0.0_pReal, & !< previous time interval
remainingLoadCaseTime = 0.0_pReal !< remaining time of current load case
logical :: &
guess, & !< guess along former trajectory
stagIterate
2019-05-15 02:14:38 +05:30
integer :: &
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
inc, & !< current increment in current load case
2019-05-15 02:14:38 +05:30
totalIncsCounter = 0, & !< total # of increments
convergedCounter = 0, & !< # of converged increments
notConvergedCounter = 0, & !< # of non-converged 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
lastRestartWritten = 0, & !< total increment # at which last restart information was written
stagIter
character(len=6) :: loadcase_string
character(len=1024) :: &
incInfo
2019-09-21 05:22:55 +05:30
type(rotation) :: R
type(tLoadCase), allocatable, dimension(:) :: loadCases !< array of all load cases
2018-08-31 12:41:09 +05:30
type(tLoadCase) :: newLoadCase
type(tSolutionState), allocatable, dimension(:) :: solres
integer(MPI_OFFSET_KIND) :: fileOffset
integer(MPI_OFFSET_KIND), dimension(:), allocatable :: outputSize
2019-05-15 02:14:38 +05:30
integer, parameter :: maxByteOut = 2147483647-4096 !< limit of one file output write https://trac.mpich.org/projects/mpich/ticket/1742
integer, parameter :: maxRealOut = maxByteOut/pReal
integer(pLongInt), dimension(2) :: outputIndex
PetscErrorCode :: ierr
2019-03-12 10:36:59 +05:30
procedure(grid_mech_spectral_basic_init), pointer :: &
mech_init
2019-03-12 10:36:59 +05:30
procedure(grid_mech_spectral_basic_forward), pointer :: &
mech_forward
2019-03-12 10:36:59 +05:30
procedure(grid_mech_spectral_basic_solution), pointer :: &
mech_solution
external :: &
quit
!--------------------------------------------------------------------------------------------------
! init DAMASK (all modules)
call CPFEM_initAll
2016-07-29 20:04:51 +05:30
write(6,'(/,a)') ' <<<+- DAMASK_spectral init -+>>>'
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
!--------------------------------------------------------------------------------------------------
! 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))
!--------------------------------------------------------------------------------------------------
! 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-03-12 10:36:59 +05:30
mech_init => grid_mech_spectral_basic_init
mech_forward => grid_mech_spectral_basic_forward
mech_solution => grid_mech_spectral_basic_solution
2019-03-23 11:22:26 +05:30
case ('polarisation')
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-03-12 23:26:28 +05:30
mech_init => grid_mech_spectral_polarisation_init
mech_forward => grid_mech_spectral_polarisation_forward
mech_solution => grid_mech_spectral_polarisation_solution
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-03-23 11:22:26 +05:30
mech_init => grid_mech_FEM_init
mech_forward => grid_mech_FEM_forward
mech_solution => grid_mech_FEM_solution
case default
2019-05-15 02:14:38 +05:30
call IO_error(error_ID = 891, ext_msg = config_numerics%getString('spectral_solver'))
end select
!--------------------------------------------------------------------------------------------------
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
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))
do
read(fileUnit, '(A)', iostat=myStat) line
2019-05-15 02:14:38 +05:30
if ( myStat /= 0) exit
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
chunkPos = IO_stringPos(line)
2019-05-15 02:14:38 +05:30
do i = 1, chunkPos(1) ! reading compulsory parameters for loadcase
select case (IO_lc(IO_stringValue(line,chunkPos,i)))
case('l','velocitygrad','velgrad','velocitygradient','fdot','dotf','f')
2019-05-15 02:14:38 +05:30
N_def = N_def + 1
case('t','time','delta')
2019-05-15 02:14:38 +05:30
N_t = N_t + 1
case('n','incs','increments','steps','logincs','logincrements','logsteps')
2019-05-15 02:14:38 +05:30
N_n = N_n + 1
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'
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
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
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
2019-05-15 02:14:38 +05:30
readIn: do i = 1, chunkPos(1)
select case (IO_lc(IO_stringValue(line,chunkPos,i)))
case('fdot','dotf','l','velocitygrad','velgrad','velocitygradient','f') ! assign values for the deformation BC matrix
temp_valueVector = 0.0_pReal
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'
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'
endif
2019-05-15 02:14:38 +05:30
do j = 1, 9
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
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
2013-01-24 00:03:46 +05:30
case('p','pk1','piolakirchhoff','stress', 's')
temp_valueVector = 0.0_pReal
2019-05-15 02:14:38 +05:30
do j = 1, 9
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
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)
case('t','time','delta') ! increment time
2019-05-15 02:14:38 +05:30
newLoadCase%time = IO_floatValue(line,chunkPos,i+1)
case('n','incs','increments','steps') ! number of increments
2019-05-15 02:14:38 +05:30
newLoadCase%incs = IO_intValue(line,chunkPos,i+1)
case('logincs','logincrements','logsteps') ! 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
case('freq','frequency','outputfreq') ! frequency of result writings
2019-05-15 02:14:38 +05:30
newLoadCase%outputfrequency = IO_intValue(line,chunkPos,i+1)
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)
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
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)))
case('deg','degree')
case('rad','radian') ! don't convert from degree to radian
2019-05-15 02:14:38 +05:30
l = 0
case default
2019-05-15 02:14:38 +05:30
k = 0
end select
2019-05-15 02:14:38 +05:30
do j = 1, 3
temp_valueVector(j) = IO_floatValue(line,chunkPos,i+k+j)
enddo
2019-09-21 05:48:09 +05:30
call R%fromEulers(temp_valueVector(1:3),degrees=(l==1))
newLoadCase%rotation = R%asMatrix()
2018-08-31 13:09:33 +05:30
case('rotation','rot') ! assign values for the rotation matrix
temp_valueVector = 0.0_pReal
2019-05-15 02:14:38 +05:30
do j = 1, 9
temp_valueVector(j) = IO_floatValue(line,chunkPos,i+j)
enddo
2019-03-08 02:50:00 +05:30
newLoadCase%rotation = math_9to33(temp_valueVector)
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
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
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
enddo
write(6,'(2x,a)') 'velocity gradient:'
2018-08-31 13:09:33 +05:30
else if (newLoadCase%deformation%myType == 'f') then
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)
else
write(6,'(2x,12a)',advance='no') ' * '
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
2018-08-31 13:09:33 +05:30
if (any(newLoadCase%stress%maskLogical .and. &
transpose(newLoadCase%stress%maskLogical) .and. &
reshape([ .false.,.true.,.true.,.true.,.false.,.true.,.true.,.true.,.false.],[ 3,3]))) &
2019-05-15 02:14:38 +05:30
errorID = 838 ! no rotation is allowed by stress BC
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
else
write(6,'(2x,12a)',advance='no') ' * '
endif
enddo; write(6,'(/)',advance='no')
enddo
if (any(abs(matmul(newLoadCase%rotation, &
2018-08-31 13:09:33 +05:30
transpose(newLoadCase%rotation))-math_I3) > &
reshape(spread(tol_math_check,1,9),[ 3,3]))&
2018-08-31 13:09:33 +05:30
.or. abs(math_det33(newLoadCase%rotation)) > &
2019-05-15 02:14:38 +05:30
1.0_pReal + tol_math_check) errorID = 846 ! given rotation matrix contains strain
2018-08-31 13:09:33 +05:30
if (any(dNeq(newLoadCase%rotation, math_I3))) &
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)
2019-05-15 02:14:38 +05:30
if (newLoadCase%time < 0.0_pReal) errorID = 834 ! negative time increment
2018-08-31 13:09:33 +05:30
write(6,'(2x,a,f12.6)') 'time: ', newLoadCase%time
2019-05-15 02:14:38 +05:30
if (newLoadCase%incs < 1) errorID = 835 ! non-positive incs count
write(6,'(2x,a,i5)') 'increments: ', newLoadCase%incs
2019-05-15 02:14:38 +05:30
if (newLoadCase%outputfrequency < 1) errorID = 836 ! non-positive result frequency
write(6,'(2x,a,i5)') '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)) &
write(6,'(2x,a,i5)') '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)
call results_openJobFile
call HDF5_closeGroup(results_addGroup('geometry'))
call results_addAttribute('grid',grid,'geometry')
call results_addAttribute('size',geomSize,'geometry')
call results_closeJobFile
!--------------------------------------------------------------------------------------------------
2018-08-31 13:09:33 +05:30
! doing initialization depending on active solvers
call Utilities_init()
do field = 1, nActiveFields
select case (loadCases(1)%ID(field))
case(FIELD_MECH_ID)
call mech_init
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
case(FIELD_DAMAGE_ID)
2019-03-12 04:21:04 +05:30
call grid_damage_spectral_init
end select
enddo
!--------------------------------------------------------------------------------------------------
! write header of output file
if (worldrank == 0) then
2019-05-15 02:14:38 +05:30
writeHeader: if (interface_restartInc < 1) then
2018-08-31 13:44:33 +05:30
open(newunit=fileUnit,file=trim(getSolverJobName())//&
'.spectralOut',form='UNFORMATTED',status='REPLACE')
write(fileUnit) 'load:', trim(loadCaseFile) ! ... and write header
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
2018-08-31 13:44:33 +05:30
write(fileUnit) 'logscales:', loadCases%logscale
write(fileUnit) 'increments:', loadCases%incs ! one entry per LoadCase
write(fileUnit) 'startingIncrement:', interface_restartInc ! start with writing out the previous inc
2018-08-31 13:44:33 +05:30
write(fileUnit) 'eoh'
close(fileUnit) ! end of header
open(newunit=statUnit,file=trim(getSolverJobName())//&
'.sta',form='FORMATTED',status='REPLACE')
write(statUnit,'(a)') 'Increment Time CutbackLevel Converged IterationsNeeded' ! statistics file
if (iand(debug_level(debug_spectral),debug_levelBasic) /= 0) &
write(6,'(/,a)') ' header of result and statistics file written out'
flush(6)
else writeHeader
open(newunit=statUnit,file=trim(getSolverJobName())//&
'.sta',form='FORMATTED', position='APPEND', status='OLD')
endif writeHeader
endif
!--------------------------------------------------------------------------------------------------
! 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
2019-05-15 02:14:38 +05:30
if (ierr /= 0) call IO_error(error_ID=894, ext_msg='MPI_allreduce')
call MPI_file_open(PETSC_COMM_WORLD, trim(getSolverJobName())//'.spectralOut', &
MPI_MODE_WRONLY + MPI_MODE_APPEND, &
MPI_INFO_NULL, &
2018-08-31 13:44:33 +05:30
fileUnit, &
ierr)
2019-05-15 02:14:38 +05:30
if (ierr /= 0) call IO_error(error_ID=894, ext_msg='MPI_file_open')
call MPI_file_get_position(fileUnit,fileOffset,ierr) ! get offset from header
if (ierr /= 0) call IO_error(error_ID=894, ext_msg='MPI_file_get_position')
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)
2019-05-15 02:14:38 +05:30
if (ierr /= 0) call IO_error(error_ID=894, ext_msg='MPI_file_seek')
2019-05-15 02:14:38 +05:30
writeUndeformed: if (interface_restartInc < 1) then
write(6,'(1/,a)') ' ... writing initial configuration to file ........................'
2019-05-15 02:14:38 +05:30
call CPFEM_results(0,0.0_pReal)
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
2019-05-15 02:14:38 +05:30
outputIndex = int([(i-1)*((maxRealOut)/materialpoint_sizeResults)+1, &
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)), &
[(outputIndex(2)-outputIndex(1)+1)*int(materialpoint_sizeResults,pLongInt)]), &
int((outputIndex(2)-outputIndex(1)+1)*int(materialpoint_sizeResults,pLongInt)), &
MPI_DOUBLE, MPI_STATUS_IGNORE, ierr)
2019-05-15 02:14:38 +05:30
if (ierr /= 0) call IO_error(error_ID=894, ext_msg='MPI_file_write')
enddo
fileOffset = fileOffset + sum(outputSize) ! forward to current file position
endif writeUndeformed
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
guess = loadCases(currentLoadCase)%followFormerTrajectory ! change of load case? homogeneous guess for the first inc
2019-05-15 02:14:38 +05:30
incLooping: do inc = 1, loadCases(currentLoadCase)%incs
totalIncsCounter = totalIncsCounter + 1
!--------------------------------------------------------------------------------------------------
! forwarding time
timeIncOld = timeinc ! last timeinc that brought former inc to an end
2019-05-15 02:14:38 +05:30
if (loadCases(currentLoadCase)%logscale == 0) then ! linear scale
timeinc = loadCases(currentLoadCase)%time/real(loadCases(currentLoadCase)%incs,pReal)
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))
endif
2018-08-31 13:09:33 +05:30
else ! not-1st load case of logarithmic scale
timeinc = time0 * &
( (1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc ,pReal)/&
real(loadCases(currentLoadCase)%incs ,pReal))&
2019-05-15 02:14:38 +05:30
-(1.0_pReal + loadCases(currentLoadCase)%time/time0 )**(real( inc-1 ,pReal)/&
real(loadCases(currentLoadCase)%incs ,pReal)))
endif
endif
timeinc = timeinc * real(subStepFactor,pReal)**real(-cutBackLevel,pReal) ! depending on cut back level, decrease time step
skipping: if (totalIncsCounter <= interface_restartInc) then ! not yet at restart inc?
time = time + timeinc ! just advance time, skip already performed calculation
guess = .true. ! QUESTION:why forced guessing instead of inheriting loadcase preference
else skipping
2019-05-15 02:14:38 +05:30
stepFraction = 0 ! fraction scaled by stepFactor**cutLevel
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
!--------------------------------------------------------------------------------------------------
! report begin of new step
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)
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
flush(6)
!--------------------------------------------------------------------------------------------------
! forward fields
do field = 1, nActiveFields
select case(loadCases(currentLoadCase)%ID(field))
case(FIELD_MECH_ID)
call mech_forward (&
guess,timeinc,timeIncOld,remainingLoadCaseTime, &
deformation_BC = loadCases(currentLoadCase)%deformation, &
stress_BC = loadCases(currentLoadCase)%stress, &
rotation_BC = loadCases(currentLoadCase)%rotation)
2019-03-12 04:07:06 +05:30
case(FIELD_THERMAL_ID); call grid_thermal_spectral_forward
2019-03-12 04:21:04 +05:30
case(FIELD_DAMAGE_ID); call grid_damage_spectral_forward
end select
enddo
!--------------------------------------------------------------------------------------------------
! solve fields
2019-05-15 02:14:38 +05:30
stagIter = 0
stagIterate = .true.
do while (stagIterate)
do field = 1, nActiveFields
select case(loadCases(currentLoadCase)%ID(field))
case(FIELD_MECH_ID)
solres(field) = mech_solution (&
incInfo,timeinc,timeIncOld, &
stress_BC = loadCases(currentLoadCase)%stress, &
rotation_BC = loadCases(currentLoadCase)%rotation)
case(FIELD_THERMAL_ID)
2019-03-12 04:07:06 +05:30
solres(field) = grid_thermal_spectral_solution(timeinc,timeIncOld,remainingLoadCaseTime)
case(FIELD_DAMAGE_ID)
2019-03-12 04:21:04 +05:30
solres(field) = grid_damage_spectral_solution(timeinc,timeIncOld,remainingLoadCaseTime)
end select
if (.not. solres(field)%converged) exit ! no solution found
enddo
2019-05-15 02:14:38 +05:30
stagIter = stagIter + 1
stagIterate = stagIter < stagItMax &
.and. all(solres(:)%converged) &
.and. .not. all(solres(:)%stagConverged) ! stationary with respect to staggered iteration
enddo
!--------------------------------------------------------------------------------------------------
! 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.
guess = .true. ! start guessing after first converged (sub)inc
if (worldrank == 0) then
write(statUnit,*) totalIncsCounter, time, cutBackLevel, &
solres%converged, solres%iterationsNeeded
flush(statUnit)
endif
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
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-08-31 13:44:33 +05:30
call MPI_file_close(fileUnit,ierr)
close(statUnit)
2019-05-15 02:14:38 +05:30
call quit(-1*(lastRestartWritten+1)) ! quit and provide information about last restart inc written
endif
enddo subStepLooping
2019-05-15 02:14:38 +05:30
cutBackLevel = max(0, cutBackLevel - 1) ! try half number of subincs next inc
if (all(solres(:)%converged)) then
2019-05-15 02:14:38 +05:30
convergedCounter = convergedCounter + 1
write(6,'(/,a,'//IO_intOut(totalIncsCounter)//',a)') & ! report converged inc
' increment ', totalIncsCounter, ' converged'
else
2019-05-15 02:14:38 +05:30
notConvergedCounter = notConvergedCounter + 1
write(6,'(/,a,'//IO_intOut(totalIncsCounter)//',a)') & ! report non-converged inc
' increment ', totalIncsCounter, ' NOT converged'
endif; flush(6)
2019-05-15 02:14:38 +05:30
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
write(6,'(1/,a)') ' ... writing results to file ......................................'
flush(6)
call materialpoint_postResults()
2018-08-31 13:44:33 +05:30
call MPI_file_seek (fileUnit,fileOffset,MPI_SEEK_SET,ierr)
2019-05-15 02:14:38 +05:30
if (ierr /= 0) call IO_error(894, ext_msg='MPI_file_seek')
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
2019-05-15 02:14:38 +05:30
outputIndex=int([(i-1)*((maxRealOut)/materialpoint_sizeResults)+1, &
2016-06-27 21:16:34 +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)),&
[(outputIndex(2)-outputIndex(1)+1)*int(materialpoint_sizeResults,pLongInt)]), &
int((outputIndex(2)-outputIndex(1)+1)*int(materialpoint_sizeResults,pLongInt)),&
MPI_DOUBLE, MPI_STATUS_IGNORE, ierr)
2019-05-15 02:14:38 +05:30
if(ierr /=0) call IO_error(894, ext_msg='MPI_file_write')
enddo
fileOffset = fileOffset + sum(outputSize) ! forward to current file position
call CPFEM_results(totalIncsCounter,time)
endif
2019-06-29 05:39:27 +05:30
if (mod(inc,loadCases(currentLoadCase)%restartFrequency) == 0) then ! at frequency of writing restart information
restartWrite = .true. ! set restart parameter for FEsolving
lastRestartWritten = inc ! QUESTION: first call to CPFEM_general will write?
endif
endif skipping
enddo incLooping
enddo loadCaseLooping
2017-07-27 19:51:02 +05:30
!--------------------------------------------------------------------------------------------------
! report summary of whole calculation
write(6,'(/,a)') ' ###########################################################################'
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)*100.0_pReal, ' %) increments converged!'
flush(6)
2018-08-31 13:44:33 +05:30
call MPI_file_close(fileUnit,ierr)
close(statUnit)
2019-05-15 02:14:38 +05:30
if (notConvergedCounter > 0) call quit(2) ! error if some are not converged
call quit(0) ! no complains ;)
end program DAMASK_spectral