Merge branch '136-skip-plasticity-calculations-for-purely-elastic' into 'development'
reduce computation effort for plasticity==none Closes #136 See merge request damask/DAMASK!452
This commit is contained in:
commit
58410709b9
|
@ -262,15 +262,16 @@ module subroutine plastic_LpAndItsTangents(Lp, dLp_dS, dLp_dFi, &
|
||||||
i, j
|
i, j
|
||||||
|
|
||||||
|
|
||||||
|
if (phase_plasticity(ph) == PLASTICITY_NONE_ID) then
|
||||||
|
Lp = 0.0_pReal
|
||||||
|
dLp_dFi = 0.0_pReal
|
||||||
|
dLp_dS = 0.0_pReal
|
||||||
|
else
|
||||||
|
|
||||||
Mp = matmul(matmul(transpose(Fi),Fi),S)
|
Mp = matmul(matmul(transpose(Fi),Fi),S)
|
||||||
|
|
||||||
|
|
||||||
plasticType: select case (phase_plasticity(ph))
|
plasticType: select case (phase_plasticity(ph))
|
||||||
|
|
||||||
case (PLASTICITY_NONE_ID) plasticType
|
|
||||||
Lp = 0.0_pReal
|
|
||||||
dLp_dMp = 0.0_pReal
|
|
||||||
|
|
||||||
case (PLASTICITY_ISOTROPIC_ID) plasticType
|
case (PLASTICITY_ISOTROPIC_ID) plasticType
|
||||||
call isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
|
call isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
|
||||||
|
|
||||||
|
@ -297,6 +298,8 @@ module subroutine plastic_LpAndItsTangents(Lp, dLp_dS, dLp_dFi, &
|
||||||
dLp_dS(i,j,1:3,1:3) = matmul(matmul(transpose(Fi),Fi),dLp_dMp(i,j,1:3,1:3)) ! ToDo: @PS: why not: dLp_dMp:(FiT Fi)
|
dLp_dS(i,j,1:3,1:3) = matmul(matmul(transpose(Fi),Fi),dLp_dMp(i,j,1:3,1:3)) ! ToDo: @PS: why not: dLp_dMp:(FiT Fi)
|
||||||
end do; end do
|
end do; end do
|
||||||
|
|
||||||
|
end if
|
||||||
|
|
||||||
end subroutine plastic_LpAndItsTangents
|
end subroutine plastic_LpAndItsTangents
|
||||||
|
|
||||||
|
|
||||||
|
@ -318,6 +321,7 @@ module function plastic_dotState(subdt,co,ip,el,ph,en) result(broken)
|
||||||
logical :: broken
|
logical :: broken
|
||||||
|
|
||||||
|
|
||||||
|
if (phase_plasticity(ph) /= PLASTICITY_NONE_ID) then
|
||||||
Mp = matmul(matmul(transpose(phase_mechanical_Fi(ph)%data(1:3,1:3,en)),&
|
Mp = matmul(matmul(transpose(phase_mechanical_Fi(ph)%data(1:3,1:3,en)),&
|
||||||
phase_mechanical_Fi(ph)%data(1:3,1:3,en)),phase_mechanical_S(ph)%data(1:3,1:3,en))
|
phase_mechanical_Fi(ph)%data(1:3,1:3,en)),phase_mechanical_S(ph)%data(1:3,1:3,en))
|
||||||
|
|
||||||
|
@ -341,8 +345,9 @@ module function plastic_dotState(subdt,co,ip,el,ph,en) result(broken)
|
||||||
case (PLASTICITY_NONLOCAL_ID) plasticType
|
case (PLASTICITY_NONLOCAL_ID) plasticType
|
||||||
call nonlocal_dotState(Mp,thermal_T(ph,en),subdt,ph,en,ip,el)
|
call nonlocal_dotState(Mp,thermal_T(ph,en),subdt,ph,en,ip,el)
|
||||||
end select plasticType
|
end select plasticType
|
||||||
broken = any(IEEE_is_NaN(plasticState(ph)%dotState(:,en)))
|
end if
|
||||||
|
|
||||||
|
broken = any(IEEE_is_NaN(plasticState(ph)%dotState(:,en)))
|
||||||
|
|
||||||
end function plastic_dotState
|
end function plastic_dotState
|
||||||
|
|
||||||
|
@ -390,8 +395,7 @@ module function plastic_deltaState(ph, en) result(broken)
|
||||||
integer, intent(in) :: &
|
integer, intent(in) :: &
|
||||||
ph, &
|
ph, &
|
||||||
en
|
en
|
||||||
logical :: &
|
logical :: broken
|
||||||
broken
|
|
||||||
|
|
||||||
real(pReal), dimension(3,3) :: &
|
real(pReal), dimension(3,3) :: &
|
||||||
Mp
|
Mp
|
||||||
|
@ -399,36 +403,35 @@ module function plastic_deltaState(ph, en) result(broken)
|
||||||
myOffset, &
|
myOffset, &
|
||||||
mySize
|
mySize
|
||||||
|
|
||||||
|
broken = .false.
|
||||||
|
|
||||||
|
select case (phase_plasticity(ph))
|
||||||
|
case (PLASTICITY_NONLOCAL_ID,PLASTICITY_KINEHARDENING_ID)
|
||||||
|
|
||||||
Mp = matmul(matmul(transpose(phase_mechanical_Fi(ph)%data(1:3,1:3,en)),&
|
Mp = matmul(matmul(transpose(phase_mechanical_Fi(ph)%data(1:3,1:3,en)),&
|
||||||
phase_mechanical_Fi(ph)%data(1:3,1:3,en)),phase_mechanical_S(ph)%data(1:3,1:3,en))
|
phase_mechanical_Fi(ph)%data(1:3,1:3,en)),&
|
||||||
|
phase_mechanical_S(ph)%data(1:3,1:3,en))
|
||||||
|
|
||||||
plasticType: select case (phase_plasticity(ph))
|
plasticType: select case (phase_plasticity(ph))
|
||||||
|
|
||||||
case (PLASTICITY_KINEHARDENING_ID) plasticType
|
case (PLASTICITY_KINEHARDENING_ID) plasticType
|
||||||
call plastic_kinehardening_deltaState(Mp,ph,en)
|
call plastic_kinehardening_deltaState(Mp,ph,en)
|
||||||
broken = any(IEEE_is_NaN(plasticState(ph)%deltaState(:,en)))
|
|
||||||
|
|
||||||
case (PLASTICITY_NONLOCAL_ID) plasticType
|
case (PLASTICITY_NONLOCAL_ID) plasticType
|
||||||
call plastic_nonlocal_deltaState(Mp,ph,en)
|
call plastic_nonlocal_deltaState(Mp,ph,en)
|
||||||
broken = any(IEEE_is_NaN(plasticState(ph)%deltaState(:,en)))
|
|
||||||
|
|
||||||
case default
|
|
||||||
broken = .false.
|
|
||||||
|
|
||||||
end select plasticType
|
end select plasticType
|
||||||
|
|
||||||
|
broken = any(IEEE_is_NaN(plasticState(ph)%deltaState(:,en)))
|
||||||
if (.not. broken) then
|
if (.not. broken) then
|
||||||
select case(phase_plasticity(ph))
|
|
||||||
case (PLASTICITY_NONLOCAL_ID,PLASTICITY_KINEHARDENING_ID)
|
|
||||||
|
|
||||||
myOffset = plasticState(ph)%offsetDeltaState
|
myOffset = plasticState(ph)%offsetDeltaState
|
||||||
mySize = plasticState(ph)%sizeDeltaState
|
mySize = plasticState(ph)%sizeDeltaState
|
||||||
plasticState(ph)%state(myOffset + 1:myOffset + mySize,en) = &
|
plasticState(ph)%state(myOffset + 1:myOffset + mySize,en) = &
|
||||||
plasticState(ph)%state(myOffset + 1:myOffset + mySize,en) + plasticState(ph)%deltaState(1:mySize,en)
|
plasticState(ph)%state(myOffset + 1:myOffset + mySize,en) + plasticState(ph)%deltaState(1:mySize,en)
|
||||||
end select
|
|
||||||
end if
|
end if
|
||||||
|
|
||||||
|
end select
|
||||||
|
|
||||||
end function plastic_deltaState
|
end function plastic_deltaState
|
||||||
|
|
||||||
|
|
||||||
|
@ -453,7 +456,7 @@ function plastic_active(plastic_label) result(active_plastic)
|
||||||
phase => phases%get(ph)
|
phase => phases%get(ph)
|
||||||
mech => phase%get('mechanical')
|
mech => phase%get('mechanical')
|
||||||
pl => mech%get('plastic',defaultVal = emptyDict)
|
pl => mech%get('plastic',defaultVal = emptyDict)
|
||||||
if(pl%get_asString('type',defaultVal='none') == plastic_label) active_plastic(ph) = .true.
|
active_plastic(ph) = pl%get_asString('type',defaultVal='none') == plastic_label
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
end function plastic_active
|
end function plastic_active
|
||||||
|
|
Loading…
Reference in New Issue