changed isotropic ductile damage to be inline with Fe..Fd..Fp decomposition

This commit is contained in:
Luv Sharma 2014-10-27 14:15:25 +00:00
parent bec0af7b06
commit b65ccb0181
3 changed files with 86 additions and 3 deletions

View File

@ -32,6 +32,7 @@ module constitutive
constitutive_getLocalDamage, &
constitutive_putLocalDamage, &
constitutive_getDamage, &
constitutive_getDamageStrain, &
constitutive_getDamageDiffusion33, &
constitutive_getAdiabaticTemperature, &
constitutive_putAdiabaticTemperature, &
@ -1186,6 +1187,62 @@ function constitutive_getTemperature(ipc, ip, el)
end select
end function constitutive_getTemperature
!--------------------------------------------------------------------------------------------------
!> @brief returns damage deformation gradient
!--------------------------------------------------------------------------------------------------
function constitutive_getDamageStrain(ipc, ip, el)
use prec, only: &
pReal
use math, only: &
math_I3
use material, only: &
material_phase, &
LOCAL_DAMAGE_none_ID, &
LOCAL_DAMAGE_brittle_ID, &
LOCAL_DAMAGE_ductile_ID, &
LOCAL_DAMAGE_gurson_ID, &
LOCAL_DAMAGE_anisotropic_ID, &
phase_damage
! use damage_brittle_iso, only: &
! constitutive_brittle_iso_getDamageStrain
! use damage_brittle_iso, only: &
! constitutive_brittle_iso_getDamageStrain
use damage_ductile, only: &
constitutive_ductile_getDamageStrain
! use damage_gurson, only: &
! constitutive_gurson_getDamageStrain
! use damage_anisotropic, only: &
! constitutive_anisotropic_getDamageStrain
implicit none
integer(pInt), intent(in) :: &
ipc, & !< grain number
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
constitutive_getDamageStrain
select case (phase_damage(material_phase(ipc,ip,el)))
case (LOCAL_DAMAGE_none_ID)
constitutive_getDamageStrain = math_I3 !
case (LOCAL_DAMAGE_brittle_ID)
constitutive_getDamageStrain = math_I3 !
case (LOCAL_DAMAGE_ductile_ID)
constitutive_getDamageStrain = constitutive_ductile_getDamageStrain(ipc, ip, el)
case (LOCAL_DAMAGE_gurson_ID)
constitutive_getDamageStrain = math_I3 !
case (LOCAL_DAMAGE_anisotropic_ID)
constitutive_getDamageStrain = math_I3 !
end select
end function constitutive_getDamageStrain
!--------------------------------------------------------------------------------------------------
!> @brief returns thermal deformation gradient

View File

@ -3297,7 +3297,8 @@ logical function crystallite_integrateStress(&
debug_StressLoopDistribution
use constitutive, only: constitutive_LpAndItsTangent, &
constitutive_TandItsTangent, &
constitutive_getThermalStrain
constitutive_getThermalStrain, &
constitutive_getDamageStrain
use math, only: math_mul33x33, &
math_mul33xx33, &
math_mul66x6, &
@ -3424,7 +3425,8 @@ logical function crystallite_integrateStress(&
endif
A = math_mul33x33(Fg_new,invFp_current) ! intermediate tensor needed later to calculate dFe_dLp
Fi = math_I3 ! intermediate configuration, assume decomposition as F = Fe Fi Fp
Fi = math_mul33x33(Fi,constitutive_getThermalStrain(g,i,e))
Fi = math_mul33x33(math_mul33x33(Fi,constitutive_getThermalStrain(g,i,e)), &
constitutive_getDamageStrain(g,i,e)) ! Fi = Ft Fd
Ci = math_mul33x33(math_transpose33(Fi),Fi) ! non-plastic stretch tensor (neglecting elastic contributions)
invFi = math_inv33(Fi)
detFi = math_det33(Fi)

View File

@ -45,6 +45,7 @@ module damage_ductile
damage_ductile_dotState, &
damage_ductile_microstructure, &
constitutive_ductile_getDamage, &
constitutive_ductile_getDamageStrain, &
constitutive_ductile_putDamage, &
damage_ductile_postResults
@ -321,7 +322,7 @@ subroutine damage_ductile_microstructure(nSlip,accumulatedSlip,ipc, ip, el)
end subroutine damage_ductile_microstructure
!--------------------------------------------------------------------------------------------------
!> @brief returns temperature based on local damage model state layout
!> @brief returns damage value based on local damage model state layout
!--------------------------------------------------------------------------------------------------
function constitutive_ductile_getDamage(ipc, ip, el)
use material, only: &
@ -339,6 +340,29 @@ function constitutive_ductile_getDamage(ipc, ip, el)
damageState(mappingConstitutive(2,ipc,ip,el))%state(1,mappingConstitutive(1,ipc,ip,el))
end function constitutive_ductile_getDamage
!--------------------------------------------------------------------------------------------------
!> @brief returns damage deformation gradient (extra intermediate configuration) based on
!> local damage model state layout
!--------------------------------------------------------------------------------------------------
function constitutive_ductile_getDamageStrain(ipc, ip, el)
use math, only: &
math_I3
use material, only: &
mappingConstitutive, &
damageState
implicit none
integer(pInt), intent(in) :: &
ipc, & !< grain number
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
constitutive_ductile_getDamageStrain
constitutive_ductile_getDamageStrain = &
math_I3 / ( &
damageState(mappingConstitutive(2,ipc,ip,el))%state(1,mappingConstitutive(1,ipc,ip,el)) )**(1_pInt/3_pInt) ! volumetric deformation gradient due to damage
end function constitutive_ductile_getDamageStrain
!--------------------------------------------------------------------------------------------------
!> @brief returns damage value based on local damage