Introducing the capability to restart jobs that crashed in the middle of sims. At the moment, this feature is exclusive for Marc.

Major changes:

CPFEM.f90 => 
1. Moving the initialization out of CPFEM_general into a separate subroutine, which is directly called by the hypela2 (Beware, the Abaqus version must also be modified in order to adapt with this change).
2. Restore primary state variables in CPFEM_init from binary files when requested (Marc flag: restart read).
3. Writing primary state variables into binary files (Marc flag: restart write).

FEsolving.f90 =>
1. Adding functions to recognize Marc restart flags: read and write and the corresponding restart file (parent job).
2. Change the initial value of cycleCounter = -1 in conjuction with the change made the ping-pong scheme

homogenization_RGC.f90 =>
1. Just syntax polishing.

IO.f90 =>
1. Adding functions/subroutines to open binary files for writing the primary state variables for restart purpose.

mpie_cpfem_marc.f90

1. Modification of the general scheme for collection and calculation in order to accommodate the newly added restart feature.
This commit is contained in:
Denny Tjahjanto 2010-11-03 14:39:18 +00:00
parent be265aef37
commit 763c20b302
5 changed files with 408 additions and 110 deletions

View File

@ -22,25 +22,166 @@ logical :: CPFEM_init_done = .false
CONTAINS CONTAINS
!*********************************************************
!*** call (thread safe) all module initializations ***
!*********************************************************
subroutine CPFEM_initAll(Temperature,element,IP)
use prec, only: pReal, &
prec_init
use numerics, only: numerics_init
use debug, only: debug_init
use FEsolving, only: FE_init
use math, only: math_init
use mesh, only: mesh_init
use lattice, only: lattice_init
use material, only: material_init
use constitutive, only: constitutive_init
use crystallite, only: crystallite_init
use homogenization, only: homogenization_init
use IO, only: IO_init
use mpie_interface
implicit none
integer(pInt), intent(in) :: element, & ! FE element number
IP ! FE integration point number
real(pReal), intent(in) :: Temperature ! temperature
real(pReal) rnd
integer(pInt) i,n
! initialization step (three dimensional stress state check missing?)
if (.not. CPFEM_init_done) then
call random_number(rnd)
do i=1,int(256.0*rnd)
n = n+1_pInt ! wasting random amount of time...
enddo ! ...to break potential race in multithreading
n = n+1_pInt
if (.not. CPFEM_init_inProgress) then ! yes my thread won!
CPFEM_init_inProgress = .true.
call prec_init()
call IO_init()
call numerics_init()
call debug_init()
call math_init()
call FE_init()
call mesh_init(IP, element) ! pass on coordinates to alter calcMode of first ip
call lattice_init()
call material_init()
call constitutive_init()
call crystallite_init(Temperature) ! (have to) use temperature of first IP for whole model
call homogenization_init(Temperature)
call CPFEM_init()
call mpie_interface_init()
CPFEM_init_done = .true.
CPFEM_init_inProgress = .false.
else ! loser, loser...
do while (CPFEM_init_inProgress)
enddo
endif
endif
end subroutine
!********************************************************* !*********************************************************
!*** allocate the arrays defined in module CPFEM *** !*** allocate the arrays defined in module CPFEM ***
!*** and initialize them *** !*** and initialize them ***
!********************************************************* !*********************************************************
subroutine CPFEM_init() subroutine CPFEM_init()
use prec, only: pInt use prec, only: pInt
use debug, only: debugger
use IO, only: IO_read_jobBinaryFile
use FEsolving, only: parallelExecution, & use FEsolving, only: parallelExecution, &
symmetricSolver symmetricSolver, &
restartRead, &
restartJob
use mesh, only: mesh_NcpElems, & use mesh, only: mesh_NcpElems, &
mesh_maxNips mesh_maxNips
use material, only: homogenization_maxNgrains, &
material_phase
use constitutive, only: constitutive_state0
use crystallite, only: crystallite_F0, &
crystallite_Fp0, &
crystallite_Lp0, &
crystallite_dPdF0, &
crystallite_Tstar0_v
use homogenization, only: homogenization_sizeState, &
homogenization_state0, &
materialpoint_F, &
materialpoint_F0
implicit none implicit none
integer(pInt) i,j,k,l,m
! initialize stress and jacobian to zero ! initialize stress and jacobian to zero
allocate(CPFEM_cs(6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_cs = 0.0_pReal allocate(CPFEM_cs(6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_cs = 0.0_pReal
allocate(CPFEM_dcsdE(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_dcsdE = 0.0_pReal allocate(CPFEM_dcsdE(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_dcsdE = 0.0_pReal
allocate(CPFEM_dcsdE_knownGood(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_dcsdE_knownGood = 0.0_pReal allocate(CPFEM_dcsdE_knownGood(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_dcsdE_knownGood = 0.0_pReal
! *** restore the last converged values of each essential variable from the binary file
if (restartRead) then
if (debugger) then
!$OMP CRITICAL (write2out)
write(6,'(a)') '<<< cpfem >>> Restored state variables of last converged step from binary files'
!$OMP END CRITICAL (write2out)
endif
if (IO_read_jobBinaryFile(777,'recordedPhase',restartJob,size(material_phase))) then
read (777,rec=1) material_phase
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedF',restartJob,size(crystallite_F0))) then
read (777,rec=1) crystallite_F0
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedFp',restartJob,size(crystallite_Fp0))) then
read (777,rec=1) crystallite_Fp0
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedLp',restartJob,size(crystallite_Lp0))) then
read (777,rec=1) crystallite_Lp0
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergeddPdF',restartJob,size(crystallite_dPdF0))) then
read (777,rec=1) crystallite_dPdF0
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedTstar',restartJob,size(crystallite_Tstar0_v))) then
read (777,rec=1) crystallite_Tstar0_v
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedStateConst',restartJob)) then
m = 0_pInt
do i = 1,homogenization_maxNgrains; do j = 1,mesh_maxNips; do k = 1,mesh_NcpElems
do l = 1,size(constitutive_state0(i,j,k)%p)
m = m+1_pInt
read(777,rec=m) constitutive_state0(i,j,k)%p(l)
enddo
enddo; enddo; enddo
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergedStateHomog',restartJob)) then
m = 0_pInt
do k = 1,mesh_NcpElems; do j = 1,mesh_maxNips
do l = 1,homogenization_sizeState(j,k)
m = m+1_pInt
read(777,rec=m) homogenization_state0(j,k)%p(l)
enddo
enddo; enddo
close (777)
endif
if (IO_read_jobBinaryFile(777,'convergeddcsdE',restartJob,size(CPFEM_dcsdE))) then
read (777,rec=1) CPFEM_dcsdE
close (777)
endif
restartRead = .false.
endif
! *** end of restoring
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
write(6,*) write(6,*)
write(6,*) '<<<+- cpfem init -+>>>' write(6,*) '<<<+- cpfem init -+>>>'
@ -69,21 +210,17 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
!*** variables and functions from other modules ***! !*** variables and functions from other modules ***!
use prec, only: pReal, & use prec, only: pReal, &
pInt, & pInt
prec_init use numerics, only: relevantStrain, &
use numerics, only: numerics_init, &
relevantStrain, &
defgradTolerance, & defgradTolerance, &
iJacoStiffness iJacoStiffness
use debug, only: debug_init, & use debug, only: debug_g, &
debug_g, &
debug_i, & debug_i, &
debug_e, & debug_e, &
debugger, & debugger, &
selectiveDebugger, & selectiveDebugger, &
verboseDebugger verboseDebugger
use FEsolving, only: FE_init, & use FEsolving, only: parallelExecution, &
parallelExecution, &
outdatedFFN1, & outdatedFFN1, &
terminallyIll, & terminallyIll, &
cycleCounter, & cycleCounter, &
@ -91,28 +228,24 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
theTime, & theTime, &
theDelta, & theDelta, &
FEsolving_execElem, & FEsolving_execElem, &
FEsolving_execIP FEsolving_execIP, &
use math, only: math_init, & restartWrite
math_identity2nd, & use math, only: math_identity2nd, &
math_mul33x33, & math_mul33x33, &
math_det3x3, & math_det3x3, &
math_I3, & math_I3, &
math_Mandel3333to66, & math_Mandel3333to66, &
math_Mandel33to6 math_Mandel33to6
use mesh, only: mesh_init, & use mesh, only: mesh_FEasCP, &
mesh_FEasCP, &
mesh_NcpElems, & mesh_NcpElems, &
mesh_maxNips, & mesh_maxNips, &
mesh_element, & mesh_element, &
FE_Nips FE_Nips
use lattice, only: lattice_init use material, only: homogenization_maxNgrains, &
use material, only: material_init, & microstructure_elemhomo, &
homogenization_maxNgrains, & material_phase
microstructure_elemhomo use constitutive, only: constitutive_state0,constitutive_state
use constitutive, only: constitutive_init,& use crystallite, only: crystallite_F0, &
constitutive_state0,constitutive_state
use crystallite, only: crystallite_init, &
crystallite_F0, &
crystallite_partionedF, & crystallite_partionedF, &
crystallite_Fp0, & crystallite_Fp0, &
crystallite_Fp, & crystallite_Fp, &
@ -122,8 +255,7 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
crystallite_dPdF, & crystallite_dPdF, &
crystallite_Tstar0_v, & crystallite_Tstar0_v, &
crystallite_Tstar_v crystallite_Tstar_v
use homogenization, only: homogenization_init, & use homogenization, only: homogenization_sizeState, &
homogenization_sizeState, &
homogenization_state, & homogenization_state, &
homogenization_state0, & homogenization_state0, &
materialpoint_F, & materialpoint_F, &
@ -134,7 +266,7 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
materialpoint_Temperature, & materialpoint_Temperature, &
materialpoint_stressAndItsTangent, & materialpoint_stressAndItsTangent, &
materialpoint_postResults materialpoint_postResults
use IO, only: IO_init use IO, only: IO_write_jobBinaryFile
use mpie_interface use mpie_interface
implicit none implicit none
@ -184,38 +316,6 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
! CPFEM_odd_stress, & ! CPFEM_odd_stress, &
! CPFEM_odd_jacobian ! CPFEM_odd_jacobian
! initialization step (three dimensional stress state check missing?)
if (.not. CPFEM_init_done) then
call random_number(rnd)
do i=1,int(256.0*rnd)
n = n+1_pInt ! wasting random amount of time...
enddo ! ...to break potential race in multithreading
n = n+1_pInt
if (.not. CPFEM_init_inProgress) then ! yes my thread won!
CPFEM_init_inProgress = .true.
call prec_init()
call IO_init()
call numerics_init()
call debug_init()
call math_init()
call FE_init()
call mesh_init(IP, element) ! pass on coordinates to alter calcMode of first ip
call lattice_init()
call material_init()
call constitutive_init()
call crystallite_init(Temperature) ! (have to) use temperature of first IP for whole model
call homogenization_init(Temperature)
call CPFEM_init()
call mpie_interface_init()
CPFEM_init_done = .true.
CPFEM_init_inProgress = .false.
else ! loser, loser...
do while (CPFEM_init_inProgress)
enddo
endif
endif
cp_en = mesh_FEasCP('elem',element) cp_en = mesh_FEasCP('elem',element)
selectiveDebugger = (cp_en == debug_e .and. IP == debug_i) selectiveDebugger = (cp_en == debug_e .and. IP == debug_i)
@ -262,6 +362,65 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
homogenization_state0(j,k)%p = homogenization_state(j,k)%p ! internal state of homogenization scheme homogenization_state0(j,k)%p = homogenization_state(j,k)%p ! internal state of homogenization scheme
enddo enddo
enddo enddo
! *** dump the last converged values of each essential variable to a binary file
if (restartWrite) then
if (debugger) then
!$OMP CRITICAL (write2out)
write(6,'(a)') '<<< cpfem >>> Writing state variables of last converged step to binary files'
!$OMP END CRITICAL (write2out)
endif
if (IO_write_jobBinaryFile(777,'recordedPhase',size(material_phase))) then
write (777,rec=1) material_phase
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedF',size(crystallite_F0))) then
write (777,rec=1) crystallite_F0
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedFp',size(crystallite_Fp0))) then
write (777,rec=1) crystallite_Fp0
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedLp',size(crystallite_Lp0))) then
write (777,rec=1) crystallite_Lp0
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergeddPdF',size(crystallite_dPdF0))) then
write (777,rec=1) crystallite_dPdF0
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedTstar',size(crystallite_Tstar0_v))) then
write (777,rec=1) crystallite_Tstar0_v
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedStateConst')) then
m = 0_pInt
do i = 1,homogenization_maxNgrains; do j = 1,mesh_maxNips; do k = 1,mesh_NcpElems
do l = 1,size(constitutive_state0(i,j,k)%p)
m = m+1_pInt
write(777,rec=m) constitutive_state0(i,j,k)%p(l)
enddo
enddo; enddo; enddo
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergedStateHomog')) then
m = 0_pInt
do k = 1,mesh_NcpElems; do j = 1,mesh_maxNips
do l = 1,homogenization_sizeState(j,k)
m = m+1_pInt
write(777,rec=m) homogenization_state0(j,k)%p(l)
enddo
enddo; enddo
close (777)
endif
if (IO_write_jobBinaryFile(777,'convergeddcsdE',size(CPFEM_dcsdE))) then
write (777,rec=1) CPFEM_dcsdE
close (777)
endif
endif
! *** end of dumping
endif endif
if (mode == 8 .or. mode == 9) then ! Abaqus explicit skips collect if (mode == 8 .or. mode == 9) then ! Abaqus explicit skips collect
@ -275,6 +434,9 @@ subroutine CPFEM_general(mode, ffn, ffn1, Temperature, dt, element, IP, cauchySt
if (.not. terminallyIll .and. .not. outdatedFFN1) then if (.not. terminallyIll .and. .not. outdatedFFN1) then
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
write(6,'(a32,x,i5,x,i2)') '°°° OUTDATED at element ip',cp_en,IP write(6,'(a32,x,i5,x,i2)') '°°° OUTDATED at element ip',cp_en,IP
write(6,'(a32,/,3(3(f10.6,x),/))') ' FFN was:',materialpoint_F(:,1,IP,cp_en), &
materialpoint_F(:,2,IP,cp_en), &
materialpoint_F(:,3,IP,cp_en)
write(6,'(a32,/,3(3(f10.6,x),/))') ' FFN1 now:',ffn1(:,1),ffn1(:,2),ffn1(:,3) write(6,'(a32,/,3(3(f10.6,x),/))') ' FFN1 now:',ffn1(:,1),ffn1(:,2),ffn1(:,3)
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
outdatedFFN1 = .true. outdatedFFN1 = .true.

View File

@ -6,21 +6,24 @@
use prec, only: pInt,pReal use prec, only: pInt,pReal
implicit none implicit none
integer(pInt) cycleCounter, theInc integer(pInt) :: cycleCounter = 0_pInt, theInc = -1_pInt
real(pReal) theTime, theDelta real(pReal) :: theTime = 0.0_pReal, theDelta = 0.0_pReal
logical :: lastIncConverged = .false.,outdatedByNewInc = .false.,outdatedFFN1 = .false.,terminallyIll = .false. logical :: lastIncConverged = .false.,outdatedByNewInc = .false.,outdatedFFN1 = .false.,terminallyIll = .false.
logical :: symmetricSolver = .false. logical :: symmetricSolver = .false.
logical :: parallelExecution = .true. logical :: parallelExecution = .true.
logical :: restartWrite = .false.
logical :: restartRead = .false.
logical :: lastMode = .true., cutBack = .false. logical :: lastMode = .true., cutBack = .false.
logical, dimension(:,:), allocatable :: calcMode logical, dimension(:,:), allocatable :: calcMode
integer(pInt), dimension(:,:), allocatable :: FEsolving_execIP integer(pInt), dimension(:,:), allocatable :: FEsolving_execIP
integer(pInt), dimension(2) :: FEsolving_execElem integer(pInt), dimension(2) :: FEsolving_execElem
character(len=1024) restartJob
CONTAINS CONTAINS
!*********************************************************** !***********************************************************
! determine wether a symmetric solver is used ! determine wether a symmetric solver is used
! and whether restart is requested
!*********************************************************** !***********************************************************
subroutine FE_init() subroutine FE_init()
@ -29,27 +32,30 @@
implicit none implicit none
integer(pInt), parameter :: fileunit = 222 integer(pInt), parameter :: fileunit = 222
integer(pInt), parameter :: maxNchunks = 2 integer(pInt), parameter :: maxNchunks = 6
integer(pInt), dimension(1+2*maxNchunks) :: positions integer(pInt), dimension(1+2*maxNchunks) :: positions
character(len=64) tag
character(len=1024) line character(len=1024) line
write(6,*)
write(6,*) '<<<+- FEsolving init -+>>>'
write(6,*) '$Id$'
write(6,*)
if (IO_open_inputFile(fileunit)) then if (IO_open_inputFile(fileunit)) then
rewind(fileunit) rewind(fileunit)
do do
read (fileunit,'(a1024)',END=100) line read (fileunit,'(a1024)',END=100) line
positions(1:1+2*1) = IO_stringPos(line,1) positions = IO_stringPos(line,maxNchunks)
if( IO_lc(IO_stringValue(line,positions,1)) == 'solver' ) then tag = IO_lc(IO_stringValue(line,positions,1)) ! extract key
read (fileunit,'(a1024)',END=100) line ! Garbage line select case(tag)
positions(1:1+2*2) = IO_stringPos(line,2) case ('solver')
read (fileunit,'(a1024)',END=100) line ! next line
positions = IO_stringPos(line,maxNchunks)
symmetricSolver = (IO_intValue(line,positions,2) /= 1_pInt) symmetricSolver = (IO_intValue(line,positions,2) /= 1_pInt)
exit case ('restart')
endif read (fileunit,'(a1024)',END=100) line ! next line
positions = IO_stringPos(line,maxNchunks)
restartWrite = iand(IO_intValue(line,positions,1),1_pInt) > 0_pInt
restartRead = iand(IO_intValue(line,positions,1),2_pInt) > 0_pInt
end select
enddo enddo
else else
call IO_error(101) ! cannot open input file call IO_error(101) ! cannot open input file
@ -57,6 +63,33 @@
100 close(fileunit) 100 close(fileunit)
if (restartRead .and. IO_open_logFile(fileunit)) then
rewind(fileunit)
do
read (fileunit,'(a1024)',END=200) line
positions = IO_stringPos(line,maxNchunks)
if ( IO_lc(IO_stringValue(line,positions,1)) == 'restart' .and. &
IO_lc(IO_stringValue(line,positions,2)) == 'file' .and. &
IO_lc(IO_stringValue(line,positions,3)) == 'job' .and. &
IO_lc(IO_stringValue(line,positions,4)) == 'id' ) &
restartJob = IO_StringValue(line,positions,6)
enddo
endif
200 close(fileunit)
!$OMP CRITICAL (write2out)
write(6,*)
write(6,*) '<<<+- FEsolving init -+>>>'
write(6,*) '$Id$'
write(6,*)
write(6,*) 'restart writing: ', restartWrite
write(6,*) 'restart reading: ', restartRead
if (restartRead) write(6,*) 'restart Job: ', trim(restartJob)
write(6,*)
!$OMP END CRITICAL (write2out)
return return
end subroutine end subroutine

View File

@ -184,6 +184,28 @@ end function
endfunction endfunction
!********************************************************************
! open FEM logfile to given unit
!********************************************************************
logical function IO_open_logFile(unit)
use prec, only: pReal, pInt
use mpie_interface
implicit none
integer(pInt), intent(in) :: unit
IO_open_logFile = .false.
open(unit,status='old',err=100,file=trim(getSolverWorkingDirectoryName())//&
trim(getSolverJobName())//LogFileExtension)
IO_open_logFile = .true.
100 return
endfunction
!******************************************************************** !********************************************************************
! open (write) file related to current job ! open (write) file related to current job
! but with different extension to given unit ! but with different extension to given unit
@ -208,6 +230,68 @@ end function
endfunction endfunction
!********************************************************************
! open (write) binary file related to current job
! but with different extension to given unit
!********************************************************************
logical function IO_write_jobBinaryFile(unit,newExt,recMultiplier)
use prec, only: pReal, pInt
use mpie_interface
implicit none
integer(pInt), intent(in) :: unit
integer(pInt), intent(in), optional :: recMultiplier
character(*), intent(in) :: newExt
IO_write_jobBinaryFile = .false.
if (present(recMultiplier)) then
open(unit,status='replace',form='unformatted',access='direct',recl=pReal*recMultiplier, &
err=100,file=trim(getSolverWorkingDirectoryName())//&
trim(getSolverJobName())//'.'//newExt)
else
open(unit,status='replace',form='unformatted',access='direct',recl=pReal, &
err=100,file=trim(getSolverWorkingDirectoryName())//&
trim(getSolverJobName())//'.'//newExt)
endif
IO_write_jobBinaryFile = .true.
100 return
endfunction
!********************************************************************
! open (read) binary file related to restored job
! and with different extension to given unit
!********************************************************************
logical function IO_read_jobBinaryFile(unit,newExt,jobName,recMultiplier)
use prec, only: pReal, pInt
use mpie_interface
implicit none
integer(pInt), intent(in) :: unit
integer(pInt), intent(in), optional :: recMultiplier
character(*), intent(in) :: newExt, jobName
IO_read_jobBinaryFile = .false.
if (present(recMultiplier)) then
open(unit,status='old',form='unformatted',access='direct',recl=pReal*recMultiplier, &
err=100,file=trim(getSolverWorkingDirectoryName())//&
trim(jobName)//'.'//newExt)
else
open(unit,status='old',form='unformatted',access='direct',recl=pReal, &
err=100,file=trim(getSolverWorkingDirectoryName())//&
trim(jobName)//'.'//newExt)
endif
IO_read_jobBinaryFile = .true.
100 return
endfunction
!******************************************************************** !********************************************************************
! hybrid IA repetition counter ! hybrid IA repetition counter
!******************************************************************** !********************************************************************

View File

@ -168,12 +168,12 @@ subroutine homogenization_RGC_init(&
mySize = 1 mySize = 1
case('volumediscrepancy') case('volumediscrepancy')
mySize = 1 mySize = 1
case default
mySize = 0
case('averagerelaxrate') case('averagerelaxrate')
mySize = 1 mySize = 1
case('maximumrelaxrate') case('maximumrelaxrate')
mySize = 1 mySize = 1
case default
mySize = 0
end select end select
if (mySize > 0_pInt) then ! any meaningful output found if (mySize > 0_pInt) then ! any meaningful output found

View File

@ -43,16 +43,19 @@ MODULE mpie_interface
character(len=64), parameter :: FEsolver = 'Marc' character(len=64), parameter :: FEsolver = 'Marc'
character(len=4), parameter :: InputFileExtension = '.dat' character(len=4), parameter :: InputFileExtension = '.dat'
character(len=4), parameter :: LogFileExtension = '.log'
CONTAINS CONTAINS
subroutine mpie_interface_init() subroutine mpie_interface_init()
!$OMP CRITICAL (write2out)
write(6,*) write(6,*)
write(6,*) '<<<+- mpie_cpfem_marc init -+>>>' write(6,*) '<<<+- mpie_cpfem_marc init -+>>>'
write(6,*) '$Id$' write(6,*) '$Id$'
write(6,*) write(6,*)
write(6,*) !$OMP END CRITICAL (write2out)
write(6,*)
return return
end subroutine end subroutine
@ -198,7 +201,7 @@ subroutine hypela2(&
use debug, only: debug_info, & use debug, only: debug_info, &
debug_reset debug_reset
use mesh, only: mesh_FEasCP use mesh, only: mesh_FEasCP
use CPFEM, only: CPFEM_general,CPFEM_init_done use CPFEM, only: CPFEM_initAll,CPFEM_general,CPFEM_init_done
implicit none implicit none
! ** Start of generated type statements ** ! ** Start of generated type statements **
@ -226,46 +229,62 @@ subroutine hypela2(&
integer(pInt) computationMode, i, cp_en integer(pInt) computationMode, i, cp_en
if ( .not. CPFEM_init_done ) then if (.not. CPFEM_init_done) call CPFEM_initAll(t(1),n(1),nn)
computationMode = 2 ! calc + init if (lovl == 4) then ! Marc requires stiffness in separate call
theTime = cptim ! record current starting time
theDelta = timinc ! record current time increment
theInc = inc ! record current increment number
!$OMP CRITICAL (write2out)
write (6,'(a,x,i6,x,i2)') '<< hypela2 >> first call special case..!',n(1),nn; call flush(6)
!$OMP END CRITICAL (write2out)
else if (lovl == 4) then ! Marc requires stiffness in separate call
if ( timinc < theDelta .and. theInc == inc ) then ! first after cutback if ( timinc < theDelta .and. theInc == inc ) then ! first after cutback
computationMode = 7 ! --> restore tangent and return computationMode = 7 ! --> restore tangent and return it
else else
computationMode = 6 ! --> just return known value computationMode = 6 ! --> just return known tangent
endif endif
else else ! stress requested (lovl == 6)
cp_en = mesh_FEasCP('elem',n(1)) cp_en = mesh_FEasCP('elem',n(1))
if (theTime < cptim .or. theInc /= inc) then ! reached convergence
lastIncConverged = .true. if (cptim > theTime .or. inc /= theInc) then ! reached convergence
outdatedByNewInc = .true.
terminallyIll = .false. terminallyIll = .false.
cycleCounter = 0 cycleCounter = 0
if (inc == 0) then ! start of analysis
lastIncConverged = .false. ! no Jacobian backup
outdatedByNewInc = .false. ! no aging of state
lastMode = .false. ! pretend last step was collection
calcMode = .false. ! pretend last step was collection
!$OMP CRITICAL (write2out)
write (6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> start of analysis..!'; call flush(6)
!$OMP END CRITICAL (write2out)
else if (inc - theInc > 1) then ! restart of broken analysis
lastIncConverged = .false. ! no Jacobian backup
outdatedByNewInc = .false. ! no aging of state
lastMode = .true. ! pretend last step was calculation
calcMode = .true. ! pretend last step was calculation
!$OMP CRITICAL (write2out)
write (6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> restart of analysis..!'; call flush(6)
!$OMP END CRITICAL (write2out)
else ! just a new inc
lastIncConverged = .true. ! Jacobian backup
outdatedByNewInc = .true. ! aging of state
lastMode = .true. ! assure last step was calculation
calcMode = .true. ! assure last step was calculation
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
write (6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> former increment converged..!'; call flush(6) write (6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> former increment converged..!'; call flush(6)
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
endif
else if ( timinc < theDelta ) then ! cutBack else if ( timinc < theDelta ) then ! cutBack
calcMode = .true. ! pretend last step was calculation
terminallyIll = .false. terminallyIll = .false.
cycleCounter = 0 cycleCounter = 0
calcMode = .true. ! pretend last step was calculation
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
write(6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> cutback detected..!'; call flush(6) write(6,'(i6,x,i2,x,a)') n(1),nn,'<< hypela2 >> cutback detected..!'; call flush(6)
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
endif
endif ! convergence treatment end
calcMode(nn,cp_en) = .not. calcMode(nn,cp_en) ! ping pong (calc <--> collect) calcMode(nn,cp_en) = .not. calcMode(nn,cp_en) ! ping pong (calc <--> collect)
if ( calcMode(nn,cp_en) ) then ! now --- CALC --- if ( calcMode(nn,cp_en) ) then ! now --- CALC ---
if ( lastMode .neqv. calcMode(nn,cp_en) ) then ! first after ping pong if ( lastMode /= calcMode(nn,cp_en) ) then ! first after ping pong
call debug_reset() ! resets debugging call debug_reset() ! resets debugging
outdatedFFN1 = .false. outdatedFFN1 = .false.
cycleCounter = cycleCounter + 1 cycleCounter = cycleCounter + 1
@ -277,7 +296,7 @@ subroutine hypela2(&
computationMode = 2 ! plain calc computationMode = 2 ! plain calc
endif endif
else ! now --- COLLECT --- else ! now --- COLLECT ---
if ( lastMode .neqv. calcMode(nn,cp_en) .and. & if ( lastMode /= calcMode(nn,cp_en) .and. &
.not. terminallyIll ) call debug_info() ! first after ping pong reports (meaningful) debugging .not. terminallyIll ) call debug_info() ! first after ping pong reports (meaningful) debugging
if ( lastIncConverged ) then if ( lastIncConverged ) then
lastIncConverged = .false. lastIncConverged = .false.