removed variable lastMode from FE_solving, used to be used for detection of first call after ping pong; now this is done by checking for a change in the lovl (macro) or a change in the calculation step (abaqus)

This commit is contained in:
Christoph Kords 2013-08-02 13:28:50 +00:00
parent 74791a6686
commit 8d6b840802
4 changed files with 45 additions and 51 deletions

View File

@ -128,7 +128,6 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
cycleCounter, &
theInc, &
calcMode, &
lastMode, &
theTime, &
theDelta, &
lastIncConverged, &
@ -217,7 +216,7 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
real(pReal) :: temperature ! temp by Abaqus is intent(in)
real(pReal), dimension(6) :: stress_h
real(pReal), dimension(6,6) :: ddsdde_h
integer(pInt) :: computationMode, i, cp_en
integer(pInt) :: computationMode, i, cp_en, lastStep
logical :: cutBack
temperature = temp ! temp is intent(in)
@ -236,21 +235,19 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
if (.not. CPFEM_init_done) call CPFEM_initAll(temperature,noel,npt)
cp_en = mesh_FEasCP('elem',noel)
if (time(2) > theTime .or. kinc /= theInc) then ! reached convergence
if (time(2) > theTime .or. kInc /= theInc) then ! reached convergence
terminallyIll = .false.
cycleCounter = -1 ! first calc step increments this to cycle = 0
if (kinc == 1) then ! start of analysis
if (kInc == 1) 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,'(i8,1x,i2,1x,a)') noel,npt,'<< UMAT >> start of analysis..!'; flush(6)
!$OMP END CRITICAL (write2out)
else if (kinc - theInc > 1) then ! restart of broken analysis
else if (kInc - 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,'(i8,1x,i2,1x,a)') noel,npt,'<< UMAT >> restart of analysis..!'; flush(6)
@ -258,7 +255,6 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
else ! just the next inc
lastIncConverged = .true. ! request Jacobian backup
outdatedByNewInc = .true. ! request aging of state
lastMode = .true. ! assure last step was calculation
calcMode = .true. ! assure last step was calculation
!$OMP CRITICAL (write2out)
write (6,'(i8,1x,i2,1x,a)') noel,npt,'<< UMAT >> new increment..!'; flush(6)
@ -278,7 +274,7 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
if (calcMode(npt,cp_en)) then ! now calc
computationMode = CPFEM_CALCRESULTS
if ( lastMode .neqv. calcMode(npt,cp_en) ) then ! first after ping pong
if ( lastStep .neqv. kStep ) then ! first after ping pong
call debug_reset() ! resets debugging
outdatedFFN1 = .false.
cycleCounter = cycleCounter + 1
@ -289,7 +285,7 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
endif
else ! now collect
computationMode = CPFEM_COLLECT
if(lastMode .neqv. calcMode(npt,cp_en) .and. .not. terminallyIll) then
if(lastStep .neqv. kStep .and. .not. terminallyIll) then
call debug_info() ! first after ping pong reports debugging
endif
if (lastIncConverged) then
@ -304,8 +300,8 @@ subroutine UMAT(STRESS,STATEV,DDSDDE,SSE,SPD,SCD,&
theTime = time(2) ! record current starting time
theDelta = dtime ! record current time increment
theInc = kinc ! record current increment number
lastMode = calcMode(npt,cp_en) ! record calculationMode
theInc = kInc ! record current increment number
lastStep = kStep ! record step number
if (iand(debug_level(debug_abaqus),debug_levelBasic) /= 0) then
!$OMP CRITICAL (write2out)

View File

@ -179,7 +179,6 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
cycleCounter, &
theInc, &
calcMode, &
lastMode, &
theTime, &
theDelta, &
lastIncConverged, &
@ -274,7 +273,7 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
logical :: cutBack
real(pReal), dimension(6) :: stress
real(pReal), dimension(6,6) :: ddsdde
integer(pInt) :: computationMode, i, cp_en, node, CPnodeID
integer(pInt) :: computationMode, i, cp_en, node, CPnodeID, lastLovl
!$ integer(pInt) :: defaultNumThreadsInt !< default value set by Marc
if (iand(debug_level(debug_MARC),debug_LEVELBASIC) /= 0_pInt) then
@ -298,11 +297,11 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
!$ call omp_set_num_threads(DAMASK_NumThreadsInt) ! set number of threads for parallel execution set by DAMASK_NUM_THREADS
computationMode = 0_pInt ! save initialization value, since does not result in any calculation
computationMode = 0_pInt ! save initialization value, since it does not result in any calculation
if (lovl == 4 ) then ! jacobian requested by marc
if (timinc < theDelta .and. theInc == inc) &
if (timinc < theDelta .and. theInc == inc .and. lastLovl /= lovl) &
computationMode = CPFEM_RESTOREJACOBIAN ! first after cutback
else ! (lovl == 6) ! stress requested by marc
elseif (lovl == 6) then ! stress requested by marc
cp_en = mesh_FEasCP('elem',m(1))
if (cptim > theTime .or. inc /= theInc) then ! reached "convergence"
terminallyIll = .false.
@ -310,8 +309,8 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
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
lastLovl = lovl ! pretend that this is NOT the first after a lovl change
!$OMP CRITICAL (write2out)
write(6,'(a,i6,1x,i2)') '<< HYPELA2 >> start of analysis..! ',m(1),nn
flush(6)
@ -319,7 +318,6 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
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,'(a,i6,1x,i2)') '<< HYPELA2 >> restart of analysis..! ',m(1),nn
@ -328,7 +326,6 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
else ! >> just the next inc <<
lastIncConverged = .true. ! request Jacobian backup
outdatedByNewInc = .true. ! request aging of state
lastMode = .true. ! assure last step was calculation
calcMode = .true. ! assure last step was calculation
!$OMP CRITICAL (write2out)
write(6,'(a,i6,1x,i2)') '<< HYPELA2 >> new increment..! ',m(1),nn
@ -345,40 +342,44 @@ subroutine hypela2(d,g,e,de,s,t,dt,ngens,m,nn,kcus,matus,ndi,nshear,disp, &
!$OMP END CRITICAL (write2out)
endif ! convergence treatment end
calcMode(nn,cp_en) = .not. calcMode(nn,cp_en) ! ping pong (calc <--> collect)
if ( calcMode(nn,cp_en) ) then ! now --- CALC ---
if ( lastMode /= calcMode(nn,cp_en) ) then ! first after ping pong
call debug_reset() ! resets debugging
outdatedFFN1 = .false.
cycleCounter = cycleCounter + 1_pInt
mesh_cellnode = mesh_build_cellnodes(mesh_node,mesh_Ncellnodes) ! update cell node coordinates
call mesh_build_ipCoordinates() ! update ip coordinates
if (usePingPong) then
calcMode(nn,cp_en) = .not. calcMode(nn,cp_en) ! ping pong (calc <--> collect)
if (calcMode(nn,cp_en)) then ! now --- CALC ---
if (lastLovl /= lovl) then ! first after ping pong
call debug_reset() ! resets debugging
outdatedFFN1 = .false.
cycleCounter = cycleCounter + 1_pInt
mesh_cellnode = mesh_build_cellnodes(mesh_node,mesh_Ncellnodes) ! update cell node coordinates
call mesh_build_ipCoordinates() ! update ip coordinates
endif
if (outdatedByNewInc) then
computationMode = ior(CPFEM_CALCRESULTS,CPFEM_AGERESULTS)
outdatedByNewInc = .false. ! reset flag
else
computationMode = CPFEM_CALCRESULTS
endif
else ! now --- COLLECT ---
if (lastLovl /= lovl .and. & .not. terminallyIll) &
call debug_info() ! first after ping pong reports (meaningful) debugging
if (lastIncConverged) then
computationMode = ior(CPFEM_COLLECT,CPFEM_BACKUPJACOBIAN) ! collect and backup Jacobian after convergence
lastIncConverged = .false. ! reset flag
else
computationMode = CPFEM_COLLECT ! plain collect
endif
do node = 1,FE_Nnodes(mesh_element(2,cp_en))
CPnodeID = mesh_element(4_pInt+node,cp_en)
mesh_node(1:3,CPnodeID) = mesh_node0(1:3,CPnodeID) + numerics_unitlength * dispt(1:3,node)
enddo
endif
if ( outdatedByNewInc ) then
computationMode = ior(CPFEM_CALCRESULTS,CPFEM_AGERESULTS)
outdatedByNewInc = .false. ! reset flag
else
computationMode = CPFEM_CALCRESULTS
endif
else ! now --- COLLECT ---
if ( lastMode /= calcMode(nn,cp_en) .and. & .not. terminallyIll) call debug_info() ! first after ping pong reports (meaningful) debugging
if ( lastIncConverged ) then
computationMode = ior(CPFEM_COLLECT,CPFEM_BACKUPJACOBIAN) ! collect and backup Jacobian after convergence
lastIncConverged = .false. ! reset flag
else
computationMode = CPFEM_COLLECT ! plain collect
endif
do node = 1,FE_Nnodes(mesh_element(2,cp_en))
CPnodeID = mesh_element(4_pInt+node,cp_en)
mesh_node(1:3,CPnodeID) = mesh_node0(1:3,CPnodeID) + numerics_unitlength * dispt(1:3,node)
enddo
endif
theTime = cptim ! record current starting time
theDelta = timinc ! record current time increment
theInc = inc ! record current increment number
lastMode = calcMode(nn,cp_en) ! record calculationMode
lastLovl = lovl ! record lovl
endif
call CPFEM_general(computationMode,usePingPong,ffn,ffn1,t(1),timinc,m(1),nn,stress,ddsdde)

View File

@ -46,7 +46,6 @@ module FEsolving
restartWrite = .false., & !< write current state to enable restart
restartRead = .false., & !< restart information to continue calculation from saved state
terminallyIll = .false., & !< at least one material point is terminally ill
lastMode = .true., & !< needs description
lastIncConverged = .false., & !< needs description
outdatedByNewInc = .false. !< needs description

View File

@ -510,7 +510,6 @@ subroutine mesh_init(ip,el)
FEsolving_execElem, &
FEsolving_execIP, &
calcMode, &
lastMode, &
modelName
implicit none
@ -624,7 +623,6 @@ subroutine mesh_init(ip,el)
allocate(calcMode(mesh_maxNips,mesh_NcpElems))
calcMode = .false. ! pretend to have collected what first call is asking (F = I)
calcMode(ip,mesh_FEasCP('elem',el)) = .true. ! first ip,el needs to be already pingponged to "calc"
lastMode = .true. ! and its mode is already known...
end subroutine mesh_init