separating functionality

This commit is contained in:
Martin Diehl 2020-12-23 06:58:54 +01:00
parent 916657e2f5
commit 82eb532193
5 changed files with 47 additions and 35 deletions

View File

@ -260,7 +260,7 @@ end subroutine CPFEM_general
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine CPFEM_forward subroutine CPFEM_forward
call crystallite_forward call homogenization_forward
call constitutive_forward call constitutive_forward
end subroutine CPFEM_forward end subroutine CPFEM_forward

View File

@ -97,7 +97,7 @@ end subroutine CPFEM_restartWrite
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine CPFEM_forward subroutine CPFEM_forward
call crystallite_forward call homogenization_forward
call constitutive_forward call constitutive_forward
end subroutine CPFEM_forward end subroutine CPFEM_forward

View File

@ -32,16 +32,11 @@ module constitutive
crystallite_F0, & !< def grad at start of FE inc crystallite_F0, & !< def grad at start of FE inc
crystallite_subF, & !< def grad to be reached at end of crystallite inc crystallite_subF, & !< def grad to be reached at end of crystallite inc
crystallite_subF0, & !< def grad at start of crystallite inc crystallite_subF0, & !< def grad at start of crystallite inc
!
crystallite_Fe, & !< current "elastic" def grad (end of converged time step) crystallite_Fe, & !< current "elastic" def grad (end of converged time step)
!
crystallite_subFp0,& !< plastic def grad at start of crystallite inc crystallite_subFp0,& !< plastic def grad at start of crystallite inc
!
crystallite_subFi0,& !< intermediate def grad at start of crystallite inc crystallite_subFi0,& !< intermediate def grad at start of crystallite inc
!
crystallite_Lp0, & !< plastic velocitiy grad at start of FE inc crystallite_Lp0, & !< plastic velocitiy grad at start of FE inc
crystallite_partitionedLp0, & !< plastic velocity grad at start of homog inc crystallite_partitionedLp0, & !< plastic velocity grad at start of homog inc
!
crystallite_S0, & !< 2nd Piola-Kirchhoff stress vector at start of FE inc crystallite_S0, & !< 2nd Piola-Kirchhoff stress vector at start of FE inc
crystallite_partitionedS0 !< 2nd Piola-Kirchhoff stress vector at start of homog inc crystallite_partitionedS0 !< 2nd Piola-Kirchhoff stress vector at start of homog inc
real(pReal), dimension(:,:,:,:,:), allocatable, public :: & real(pReal), dimension(:,:,:,:,:), allocatable, public :: &
@ -170,10 +165,11 @@ module constitutive
integer, intent(in) :: ph, me integer, intent(in) :: ph, me
end subroutine constitutive_mech_windForward end subroutine constitutive_mech_windForward
module subroutine constitutive_mech_forward
end subroutine constitutive_mech_forward
! == cleaned:end =================================================================================== ! == cleaned:end ===================================================================================
module subroutine source_damage_anisoBrittle_dotState(S, ipc, ip, el) module subroutine source_damage_anisoBrittle_dotState(S, ipc, ip, el)
integer, intent(in) :: & integer, intent(in) :: &
ipc, & !< component-ID of integration point ipc, & !< component-ID of integration point
@ -381,7 +377,6 @@ module constitutive
crystallite_push33ToRef, & crystallite_push33ToRef, &
crystallite_restartWrite, & crystallite_restartWrite, &
crystallite_restartRead, & crystallite_restartRead, &
crystallite_forward, &
constitutive_initializeRestorationPoints, & constitutive_initializeRestorationPoints, &
constitutive_windForward, & constitutive_windForward, &
crystallite_restore crystallite_restore
@ -778,6 +773,12 @@ subroutine constitutive_forward
integer :: i, j integer :: i, j
crystallite_F0 = crystallite_partitionedF
crystallite_Lp0 = crystallite_Lp
crystallite_S0 = crystallite_S
call constitutive_mech_forward()
do i = 1, size(sourceState) do i = 1, size(sourceState)
do j = 1,phase_Nsources(i) do j = 1,phase_Nsources(i)
sourceState(i)%p(j)%state0 = sourceState(i)%p(j)%state sourceState(i)%p(j)%state0 = sourceState(i)%p(j)%state
@ -1605,29 +1606,4 @@ subroutine crystallite_restartRead
end subroutine crystallite_restartRead end subroutine crystallite_restartRead
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
! ToDo: Any guessing for the current states possible?
!--------------------------------------------------------------------------------------------------
subroutine crystallite_forward
integer :: i, j
crystallite_F0 = crystallite_partitionedF
crystallite_Lp0 = crystallite_Lp
crystallite_S0 = crystallite_S
do i = 1, size(plasticState)
plasticState(i)%state0 = plasticState(i)%state
constitutive_mech_Fi0(i) = constitutive_mech_Fi(i)
constitutive_mech_Fp0(i) = constitutive_mech_Fp(i)
constitutive_mech_Li0(i) = constitutive_mech_Li(i)
enddo
do i = 1,size(material_name_homogenization)
homogState (i)%state0 = homogState (i)%state
damageState (i)%state0 = damageState (i)%state
enddo
end subroutine crystallite_forward
end module constitutive end module constitutive

View File

@ -1425,5 +1425,24 @@ module subroutine constitutive_mech_windForward(ph,me)
end subroutine constitutive_mech_windForward end subroutine constitutive_mech_windForward
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
! ToDo: Any guessing for the current states possible?
!--------------------------------------------------------------------------------------------------
module subroutine constitutive_mech_forward()
integer :: ph
do ph = 1, size(plasticState)
plasticState(ph)%state0 = plasticState(ph)%state
constitutive_mech_Fi0(ph) = constitutive_mech_Fi(ph)
constitutive_mech_Fp0(ph) = constitutive_mech_Fp(ph)
constitutive_mech_Li0(ph) = constitutive_mech_Li(ph)
enddo
end subroutine constitutive_mech_forward
end submodule constitutive_mech end submodule constitutive_mech

View File

@ -112,6 +112,7 @@ module homogenization
public :: & public :: &
homogenization_init, & homogenization_init, &
materialpoint_stressAndItsTangent, & materialpoint_stressAndItsTangent, &
homogenization_forward, &
homogenization_results homogenization_results
contains contains
@ -425,4 +426,20 @@ subroutine homogenization_results
end subroutine homogenization_results end subroutine homogenization_results
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
! ToDo: Any guessing for the current states possible?
!--------------------------------------------------------------------------------------------------
subroutine homogenization_forward
integer :: ho
do ho = 1, size(material_name_homogenization)
homogState (ho)%state0 = homogState (ho)%state
damageState(ho)%state0 = damageState(ho)%state
enddo
end subroutine homogenization_forward
end module homogenization end module homogenization