DAMASK_EICMD/src/mesh/DAMASK_mesh.f90

343 lines
18 KiB
Fortran
Raw Normal View History

2018-08-17 03:44:25 +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 FEM solver
2019-10-24 01:34:50 +05:30
!> @details doing cutbacking, reporting statistics, writing
2018-08-17 03:44:25 +05:30
!> results
!--------------------------------------------------------------------------------------------------
program DAMASK_mesh
2018-09-28 10:55:32 +05:30
#include <petsc/finclude/petscsys.h>
2019-06-11 13:18:07 +05:30
use PetscDM
use prec
use DAMASK_interface
2020-09-13 13:49:38 +05:30
use parallelization
2019-06-11 13:18:07 +05:30
use IO
use math
use CPFEM2
2020-08-15 19:32:10 +05:30
use config
2020-03-20 19:38:07 +05:30
use discretization_mesh
2019-06-11 13:18:07 +05:30
use FEM_Utilities
2021-02-09 03:51:53 +05:30
use mesh_mechanical_FEM
2019-06-11 13:18:07 +05:30
implicit none
2018-08-17 03:44:25 +05:30
type :: tLoadCase
real(pReal) :: time = 0.0_pReal !< length of increment
integer :: incs = 0, & !< number of increments
outputfrequency = 1 !< frequency of result writes
logical :: followFormerTrajectory = .true. !< follow trajectory of former loadcase
integer, allocatable, dimension(:) :: faceID
type(tFieldBC), allocatable, dimension(:) :: fieldBC
end type tLoadCase
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! variables related to information from load case and geom file
2019-06-11 13:18:07 +05:30
integer, allocatable, dimension(:) :: chunkPos ! this is longer than needed for geometry parsing
integer :: &
N_def = 0 !< # of rate of deformation specifiers found in load case file
2020-08-09 10:07:50 +05:30
character(len=:), allocatable :: &
2019-06-11 13:18:07 +05:30
line
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
2019-06-11 13:18:07 +05:30
integer, parameter :: &
subStepFactor = 2 !< for each substep, divide the last time increment by 2.0
real(pReal) :: &
2020-02-04 03:29:38 +05:30
time = 0.0_pReal, & !< elapsed time
time0 = 0.0_pReal, & !< begin of interval
timeinc = 0.0_pReal, & !< current time interval
timeIncOld = 0.0_pReal, & !< previous time interval
remainingLoadCaseTime = 0.0_pReal !< remaining time of current load case
2019-06-11 13:18:07 +05:30
logical :: &
2020-02-04 03:29:38 +05:30
guess, & !< guess along former trajectory
2019-06-11 13:18:07 +05:30
stagIterate
integer :: &
2020-02-04 03:29:38 +05:30
l, &
2019-06-11 13:18:07 +05:30
i, &
errorID, &
cutBackLevel = 0, & !< cut back level \f$ t = \frac{t_{inc}}{2^l} \f$
2020-02-04 03:14:57 +05:30
stepFraction = 0, & !< fraction of current time interval
2019-06-11 13:18:07 +05:30
currentLoadcase = 0, & !< current load case
currentFace = 0, &
inc, & !< current increment in current load case
totalIncsCounter = 0, & !< total # of increments
statUnit = 0, & !< file unit for statistics output
stagIter, &
component
2020-06-17 20:17:13 +05:30
class(tNode), pointer :: &
2020-06-19 06:02:33 +05:30
num_mesh
2020-02-04 03:29:38 +05:30
character(len=pStringLen), dimension(:), allocatable :: fileContent
2020-01-26 16:28:13 +05:30
character(len=pStringLen) :: &
incInfo, &
loadcase_string
2020-06-29 20:35:11 +05:30
integer :: &
stagItMax, & !< max number of field level staggered iterations
maxCutBack !< max number of cutbacks
2020-06-24 20:15:13 +05:30
2020-02-04 03:29:38 +05:30
type(tLoadCase), allocatable, dimension(:) :: loadCases !< array of all load cases
2019-06-11 13:18:07 +05:30
type(tSolutionState), allocatable, dimension(:) :: solres
PetscInt :: faceSet, currentFaceSet, dimPlex
2019-06-11 13:18:07 +05:30
PetscErrorCode :: ierr
2020-01-21 12:16:32 +05:30
integer(kind(COMPONENT_UNDEFINED_ID)) :: ID
2019-06-11 13:18:07 +05:30
external :: &
quit
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! init DAMASK (all modules)
2019-06-11 13:18:07 +05:30
call CPFEM_initAll
2020-09-22 16:39:12 +05:30
print'(/,a)', ' <<<+- DAMASK_mesh init -+>>>'; flush(IO_STDOUT)
2020-06-17 20:17:13 +05:30
!---------------------------------------------------------------------
2020-06-17 20:17:13 +05:30
! reading field information from numerics file and do sanity checks
2020-09-13 14:47:49 +05:30
num_mesh => config_numerics%get('mesh', defaultVal=emptyDict)
2020-06-29 20:35:11 +05:30
stagItMax = num_mesh%get_asInt('maxStaggeredIter',defaultVal=10)
maxCutBack = num_mesh%get_asInt('maxCutBack',defaultVal=3)
2020-06-17 20:17:13 +05:30
if (stagItMax < 0) call IO_error(301,ext_msg='maxStaggeredIter')
if (maxCutBack < 0) call IO_error(301,ext_msg='maxCutBack')
2020-06-29 18:39:13 +05:30
! reading basic information from load case file and allocate data structure containing load cases
call DMGetDimension(geomMesh,dimPlex,ierr) !< dimension of mesh (2D or 3D)
CHKERRA(ierr)
allocate(solres(1))
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! reading basic information from load case file and allocate data structure containing load cases
2020-09-13 14:48:09 +05:30
fileContent = IO_readlines(trim(interface_loadFile))
2020-02-04 03:29:38 +05:30
do l = 1, size(fileContent)
line = fileContent(l)
2019-06-11 13:18:07 +05:30
if (IO_isBlank(line)) cycle ! skip empty lines
2019-06-11 13:18:07 +05:30
chunkPos = IO_stringPos(line)
do i = 1, chunkPos(1) ! reading compulsory parameters for loadcase
select case (IO_stringValue(line,chunkPos,i))
case('$Loadcase')
2019-06-11 13:18:07 +05:30
N_def = N_def + 1
end select
enddo ! count all identifiers to allocate memory and do sanity check
enddo
if(N_def < 1) call IO_error(error_ID = 837)
allocate(loadCases(N_def))
2019-06-11 13:18:07 +05:30
do i = 1, size(loadCases)
allocate(loadCases(i)%fieldBC(1))
loadCases(i)%fieldBC(1)%ID = FIELD_MECH_ID
2019-06-11 13:18:07 +05:30
enddo
2019-06-11 13:18:07 +05:30
do i = 1, size(loadCases)
loadCases(i)%fieldBC(1)%nComponents = dimPlex !< X, Y (, Z) displacements
allocate(loadCases(i)%fieldBC(1)%componentBC(loadCases(i)%fieldBC(1)%nComponents))
do component = 1, loadCases(i)%fieldBC(1)%nComponents
select case (component)
case (1)
loadCases(i)%fieldBC(1)%componentBC(component)%ID = COMPONENT_MECH_X_ID
case (2)
loadCases(i)%fieldBC(1)%componentBC(component)%ID = COMPONENT_MECH_Y_ID
case (3)
loadCases(i)%fieldBC(1)%componentBC(component)%ID = COMPONENT_MECH_Z_ID
2019-06-11 13:18:07 +05:30
end select
enddo
do component = 1, loadCases(i)%fieldBC(1)%nComponents
allocate(loadCases(i)%fieldBC(1)%componentBC(component)%Value(mesh_Nboundaries), source = 0.0_pReal)
allocate(loadCases(i)%fieldBC(1)%componentBC(component)%Mask (mesh_Nboundaries), source = .false.)
enddo
enddo
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! reading the load case and assign values to the allocated data structure
2020-02-04 03:29:38 +05:30
do l = 1, size(fileContent)
line = fileContent(l)
2019-06-11 13:18:07 +05:30
if (IO_isBlank(line)) cycle ! skip empty lines
2019-06-11 13:18:07 +05:30
chunkPos = IO_stringPos(line)
do i = 1, chunkPos(1)
select case (IO_stringValue(line,chunkPos,i))
!--------------------------------------------------------------------------------------------------
! loadcase information
case('$Loadcase')
2019-06-11 13:18:07 +05:30
currentLoadCase = IO_intValue(line,chunkPos,i+1)
case('Face')
2019-06-11 13:18:07 +05:30
currentFace = IO_intValue(line,chunkPos,i+1)
currentFaceSet = -1
do faceSet = 1, mesh_Nboundaries
if (mesh_boundaries(faceSet) == currentFace) currentFaceSet = faceSet
enddo
if (currentFaceSet < 0) call IO_error(error_ID = 837, ext_msg = 'invalid BC')
case('t')
2019-06-11 13:18:07 +05:30
loadCases(currentLoadCase)%time = IO_floatValue(line,chunkPos,i+1)
case('N')
2019-06-11 13:18:07 +05:30
loadCases(currentLoadCase)%incs = IO_intValue(line,chunkPos,i+1)
case('f_out')
loadCases(currentLoadCase)%outputfrequency = IO_intValue(line,chunkPos,i+1)
case('estimate_rate')
2019-06-11 13:18:07 +05:30
loadCases(currentLoadCase)%followFormerTrajectory = .false. ! do not continue to predict deformation along former trajectory
!--------------------------------------------------------------------------------------------------
! boundary condition information
case('X','Y','Z')
select case(IO_stringValue(line,chunkPos,i))
case('X')
2020-01-21 12:16:32 +05:30
ID = COMPONENT_MECH_X_ID
case('Y')
2020-01-21 12:16:32 +05:30
ID = COMPONENT_MECH_Y_ID
case('Z')
2020-01-21 12:16:32 +05:30
ID = COMPONENT_MECH_Z_ID
end select
do component = 1, loadcases(currentLoadCase)%fieldBC(1)%nComponents
if (loadCases(currentLoadCase)%fieldBC(1)%componentBC(component)%ID == ID) then
loadCases(currentLoadCase)%fieldBC(1)%componentBC(component)%Mask (currentFaceSet) = &
.true.
loadCases(currentLoadCase)%fieldBC(1)%componentBC(component)%Value(currentFaceSet) = &
IO_floatValue(line,chunkPos,i+1)
2020-01-21 12:16:32 +05:30
endif
enddo
2019-06-11 13:18:07 +05:30
end select
2020-02-04 03:29:38 +05:30
enddo
enddo
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! consistency checks and output of load case
2019-06-11 13:18:07 +05:30
loadCases(1)%followFormerTrajectory = .false. ! cannot guess along trajectory for first inc of first currentLoadCase
errorID = 0
checkLoadcases: do currentLoadCase = 1, size(loadCases)
2020-01-26 16:28:13 +05:30
write (loadcase_string, '(i0)' ) currentLoadCase
2020-09-19 11:50:29 +05:30
print'(a,i0)', ' load case: ', currentLoadCase
2020-01-26 16:28:13 +05:30
if (.not. loadCases(currentLoadCase)%followFormerTrajectory) &
2020-09-19 11:50:29 +05:30
print'(a)', ' drop guessing along trajectory'
print'(a)', ' Field '//trim(FIELD_MECH_label)
do faceSet = 1, mesh_Nboundaries
do component = 1, loadCases(currentLoadCase)%fieldBC(1)%nComponents
if (loadCases(currentLoadCase)%fieldBC(1)%componentBC(component)%Mask(faceSet)) &
print'(a,i2,a,i2,a,f12.7)', ' Face ', mesh_boundaries(faceSet), &
' Component ', component, &
' Value ', loadCases(currentLoadCase)%fieldBC(1)% &
componentBC(component)%Value(faceSet)
enddo
2020-01-26 16:28:13 +05:30
enddo
2020-09-19 11:50:29 +05:30
print'(a,f12.6)', ' time: ', loadCases(currentLoadCase)%time
2020-01-26 16:28:13 +05:30
if (loadCases(currentLoadCase)%incs < 1) errorID = 835 ! non-positive incs count
2020-09-19 11:50:29 +05:30
print'(a,i5)', ' increments: ', loadCases(currentLoadCase)%incs
2020-01-26 16:28:13 +05:30
if (loadCases(currentLoadCase)%outputfrequency < 1) errorID = 836 ! non-positive result frequency
2020-09-19 11:50:29 +05:30
print'(a,i5)', ' output frequency: ', &
2020-01-26 16:28:13 +05:30
loadCases(currentLoadCase)%outputfrequency
if (errorID > 0) call IO_error(error_ID = errorID, ext_msg = loadcase_string) ! exit with error message
2019-06-11 13:18:07 +05:30
enddo checkLoadcases
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! doing initialization depending on active solvers
2020-06-17 21:32:22 +05:30
call FEM_Utilities_init
call FEM_mechanical_init(loadCases(1)%fieldBC(1))
2019-06-11 13:18:07 +05:30
2020-02-04 03:39:46 +05:30
if (worldrank == 0) then
open(newunit=statUnit,file=trim(getSolverJobName())//'.sta',form='FORMATTED',status='REPLACE')
write(statUnit,'(a)') 'Increment Time CutbackLevel Converged IterationsNeeded' ! statistics file
endif
2019-06-11 13:18:07 +05:30
2021-05-31 23:11:22 +05:30
print'(/,a)', ' ... writing initial configuration to file ........................'
flush(IO_STDOUT)
call CPFEM_results(0,0.0_pReal)
2019-06-11 13:18:07 +05:30
loadCaseLooping: do currentLoadCase = 1, size(loadCases)
time0 = time ! load case start time
guess = loadCases(currentLoadCase)%followFormerTrajectory ! change of load case? homogeneous guess for the first inc
2019-06-11 13:18:07 +05:30
incLooping: do inc = 1, loadCases(currentLoadCase)%incs
totalIncsCounter = totalIncsCounter + 1
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! forwarding time
2019-06-11 13:18:07 +05:30
timeIncOld = timeinc ! last timeinc that brought former inc to an end
timeinc = loadCases(currentLoadCase)%time/real(loadCases(currentLoadCase)%incs,pReal)
2019-06-11 13:18:07 +05:30
timeinc = timeinc * real(subStepFactor,pReal)**real(-cutBackLevel,pReal) ! depending on cut back level, decrease time step
stepFraction = 0 ! fraction scaled by stepFactor**cutLevel
2019-10-24 01:34:50 +05:30
subStepLooping: do while (stepFraction < subStepFactor**cutBackLevel)
remainingLoadCaseTime = loadCases(currentLoadCase)%time+time0 - time
time = time + timeinc ! forward target time
stepFraction = stepFraction + 1 ! count step
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! report begin of new step
2020-09-19 11:50:29 +05:30
print'(/,a)', ' ###########################################################################'
print'(1x,a,es12.5,6(a,i0))',&
2019-10-24 01:34:50 +05:30
'Time', time, &
's: Increment ', inc, '/', loadCases(currentLoadCase)%incs,&
'-', stepFraction, '/', subStepFactor**cutBackLevel,&
' of load case ', currentLoadCase,'/',size(loadCases)
write(incInfo,'(4(a,i0))') &
2019-10-24 01:34:50 +05:30
'Increment ',totalIncsCounter,'/',sum(loadCases%incs),&
'-',stepFraction, '/', subStepFactor**cutBackLevel
2020-09-22 16:39:12 +05:30
flush(IO_STDOUT)
2018-08-17 03:44:25 +05:30
call FEM_mechanical_forward(guess,timeinc,timeIncOld,loadCases(currentLoadCase)%fieldBC(1))
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! solve fields
2019-10-24 01:34:50 +05:30
stagIter = 0
stagIterate = .true.
do while (stagIterate)
solres(1) = FEM_mechanical_solution(incInfo,timeinc,timeIncOld,loadCases(currentLoadCase)%fieldBC(1))
if(.not. solres(1)%converged) exit
2019-06-11 13:18:07 +05:30
2019-10-24 01:34:50 +05:30
stagIter = stagIter + 1
2020-06-29 20:35:11 +05:30
stagIterate = stagIter < stagItMax &
2019-10-24 01:34:50 +05:30
.and. all(solres(:)%converged) &
.and. .not. all(solres(:)%stagConverged) ! stationary with respect to staggered iteration
enddo
! check solution
cutBack = .False.
2019-10-24 01:34:50 +05:30
if(.not. all(solres(:)%converged .and. solres(:)%stagConverged)) then ! no solution found
2020-06-29 20:35:11 +05:30
if (cutBackLevel < maxCutBack) then ! do cut back
2019-10-24 01:34:50 +05:30
cutBack = .True.
stepFraction = (stepFraction - 1) * subStepFactor ! adjust to new denominator
cutBackLevel = cutBackLevel + 1
time = time - timeinc ! rewind time
timeinc = timeinc/2.0_pReal
2021-03-28 04:28:49 +05:30
print'(/,a)', ' cutting back'
else ! default behavior, exit if spectral solver does not converge
2021-03-28 04:28:49 +05:30
if (worldrank == 0) close(statUnit)
call IO_error(950)
2019-06-11 13:18:07 +05:30
endif
else
2019-10-24 01:34:50 +05:30
guess = .true. ! start guessing after first converged (sub)inc
timeIncOld = timeinc
2019-06-11 13:18:07 +05:30
endif
2020-02-04 03:39:46 +05:30
if (.not. cutBack .and. worldrank == 0) &
write(statUnit,*) totalIncsCounter, time, cutBackLevel, &
2019-10-24 01:34:50 +05:30
solres%converged, solres%iterationsNeeded ! write statistics about accepted solution
enddo subStepLooping
cutBackLevel = max(0, cutBackLevel - 1) ! try half number of subincs next inc
if (all(solres(:)%converged)) then
2020-09-19 11:50:29 +05:30
print'(/,a,i0,a)', ' increment ', totalIncsCounter, ' converged'
2019-10-24 01:34:50 +05:30
else
2020-09-19 11:50:29 +05:30
print'(/,a,i0,a)', ' increment ', totalIncsCounter, ' NOT converged'
2020-09-22 16:39:12 +05:30
endif; flush(IO_STDOUT)
2019-10-24 01:34:50 +05:30
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
2020-09-19 11:50:29 +05:30
print'(/,a)', ' ... writing results to file ......................................'
call FEM_mechanical_updateCoords
2019-10-24 01:34:50 +05:30
call CPFEM_results(totalIncsCounter,time)
endif
2019-06-11 13:18:07 +05:30
2018-08-17 03:44:25 +05:30
enddo incLooping
2019-06-11 13:18:07 +05:30
enddo loadCaseLooping
2018-08-17 03:44:25 +05:30
!--------------------------------------------------------------------------------------------------
! report summary of whole calculation
2020-09-19 11:50:29 +05:30
print'(/,a)', ' ###########################################################################'
2020-02-04 03:39:46 +05:30
if (worldrank == 0) close(statUnit)
2019-06-11 13:18:07 +05:30
call quit(0) ! no complains ;)
2018-08-17 03:44:25 +05:30
end program DAMASK_mesh