small functions with one task are better

This commit is contained in:
Martin Diehl 2019-10-24 06:16:42 +02:00
parent a90e48175c
commit 0483fc7b3f
4 changed files with 75 additions and 68 deletions

View File

@ -35,7 +35,8 @@ module CPFEM2
public :: &
CPFEM_age, &
CPFEM_initAll, &
CPFEM_results
CPFEM_results, &
CPFEM_restartWrite
contains
@ -120,15 +121,11 @@ end subroutine CPFEM_init
!--------------------------------------------------------------------------------------------------
!> @brief forwards data after successful increment
!> @brief Forward data after successful increment.
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_age(restartWrite)
subroutine CPFEM_age
logical :: restartWrite
integer :: i, ph, homog, mySource
character(len=32) :: rankStr, PlasticItem, HomogItem
integer(HID_T) :: fileHandle, groupPlastic, groupHomog
integer :: i, homog, mySource
if (iand(debug_level(debug_CPFEM), debug_levelBasic) /= 0) &
write(6,'(a)') '<< CPFEM >> aging states'
@ -153,7 +150,18 @@ subroutine CPFEM_age(restartWrite)
damageState (homog)%state0 = damageState (homog)%state
enddo
if (restartWrite) then
end subroutine CPFEM_age
!--------------------------------------------------------------------------------------------------
!> @brief Store DAMASK restart data.
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_restartWrite
integer :: ph, homog
character(len=32) :: rankStr, PlasticItem, HomogItem
integer(HID_T) :: fileHandle, groupPlastic, groupHomog
if (iand(debug_level(debug_CPFEM), debug_levelBasic) /= 0) &
write(6,'(a)') '<< CPFEM >> writing restart variables of last converged step to hdf5 file'
@ -182,17 +190,12 @@ subroutine CPFEM_age(restartWrite)
call HDF5_closeGroup(groupHomog)
call HDF5_closeFile(fileHandle)
restartWrite = .false.
endif
if (iand(debug_level(debug_CPFEM), debug_levelBasic) /= 0) &
write(6,'(a)') '<< CPFEM >> done aging states'
end subroutine CPFEM_age
end subroutine CPFEM_restartWrite
!--------------------------------------------------------------------------------------------------
!> @brief triggers writing of the results
!> @brief Trigger writing of results.
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_results(inc,time)

View File

@ -301,10 +301,11 @@ subroutine grid_mech_FEM_forward(guess,timeinc,timeinc_old,loadCaseTime,deformat
call DMDAVecGetArrayF90(mech_grid,solution_lastInc,u_lastInc,ierr); CHKERRQ(ierr)
if (cutBack) then
C_volAvg = C_volAvgLastInc ! QUESTION: where is this required?
C_volAvg = C_volAvgLastInc
else
!--------------------------------------------------------------------------------------------------
! restart information for spectral solver
call CPFEM_age ! age state and kinematics
call utilities_updateCoords(F)
if (restartWrite) then
write(6,'(/,a)') ' writing converged results for restart';flush(6)
@ -324,9 +325,9 @@ subroutine grid_mech_FEM_forward(guess,timeinc,timeinc_old,loadCaseTime,deformat
call HDF5_closeFile(fileHandle)
call CPFEM_restartWrite
restartWrite = .false.
endif
call CPFEM_age(restartWrite) ! age state and kinematics
call utilities_updateCoords(F)
C_volAvgLastInc = C_volAvg

View File

@ -279,11 +279,12 @@ subroutine grid_mech_spectral_basic_forward(guess,timeinc,timeinc_old,loadCaseTi
call DMDAVecGetArrayF90(da,solution_vec,F,ierr); CHKERRQ(ierr)
if (cutBack) then
C_volAvg = C_volAvgLastInc ! QUESTION: where is this required?
C_minMaxAvg = C_minMaxAvgLastInc ! QUESTION: where is this required?
C_volAvg = C_volAvgLastInc
C_minMaxAvg = C_minMaxAvgLastInc
else
!--------------------------------------------------------------------------------------------------
! restart information for spectral solver
call CPFEM_age ! age state and kinematics
call utilities_updateCoords(F)
if (restartWrite) then
write(6,'(/,a)') ' writing converged results for restart';flush(6)
@ -301,10 +302,10 @@ subroutine grid_mech_spectral_basic_forward(guess,timeinc,timeinc_old,loadCaseTi
call HDF5_write(fileHandle,C_minMaxAvg, 'C_minMaxAvg')
call HDF5_closeFile(fileHandle)
endif
call CPFEM_age(restartWrite) ! age state and kinematics
call utilities_updateCoords(F)
call CPFEM_restartWrite
restartWrite = .false.
endif
C_volAvgLastInc = C_volAvg
C_minMaxAvgLastInc = C_minMaxAvg

View File

@ -300,11 +300,12 @@ subroutine grid_mech_spectral_polarisation_forward(guess,timeinc,timeinc_old,loa
F_tau => FandF_tau(9:17,:,:,:)
if (cutBack) then
C_volAvg = C_volAvgLastInc ! QUESTION: where is this required?
C_minMaxAvg = C_minMaxAvgLastInc ! QUESTION: where is this required?
C_volAvg = C_volAvgLastInc
C_minMaxAvg = C_minMaxAvgLastInc
else
!--------------------------------------------------------------------------------------------------
! restart information for spectral solver
call CPFEM_age ! age state and kinematics
call utilities_updateCoords(F)
if (restartWrite) then
write(6,'(/,a)') ' writing converged results for restart';flush(6)
@ -323,10 +324,10 @@ subroutine grid_mech_spectral_polarisation_forward(guess,timeinc,timeinc_old,loa
call HDF5_write(fileHandle,C_volAvgLastInc,'C_volAvgLastInc')
call HDF5_closeFile(fileHandle)
endif
call CPFEM_age(restartWrite) ! age state and kinematics
call utilities_updateCoords(F)
call CPFEM_restartWrite
restartWrite = .false.
endif
C_volAvgLastInc = C_volAvg
C_minMaxAvgLastInc = C_minMaxAvg
@ -357,6 +358,7 @@ subroutine grid_mech_spectral_polarisation_forward(guess,timeinc,timeinc_old,loa
F_lastInc = reshape(F, [3,3,grid(1),grid(2),grid3]) ! winding F forward
F_tau_lastInc = reshape(F_tau, [3,3,grid(1),grid(2),grid3]) ! winding F_tau forward
materialpoint_F0 = reshape(F_lastInc, [3,3,1,product(grid(1:2))*grid3]) ! set starting condition for materialpoint_stressAndItsTangent
endif
!--------------------------------------------------------------------------------------------------