DAMASK_EICMD/src/phase_mechanical_plastic.f90

443 lines
16 KiB
Fortran
Raw Normal View History

2021-02-13 23:27:41 +05:30
submodule(phase:mechanical) plastic
2021-01-26 05:41:32 +05:30
interface
2021-01-26 13:09:17 +05:30
2021-01-27 05:02:44 +05:30
module function plastic_none_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_none_init
module function plastic_isotropic_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_isotropic_init
module function plastic_phenopowerlaw_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_phenopowerlaw_init
module function plastic_kinehardening_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_kinehardening_init
module function plastic_dislotwin_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_dislotwin_init
module function plastic_dislotungsten_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_dislotungsten_init
module function plastic_nonlocal_init() result(myPlasticity)
logical, dimension(:), allocatable :: &
myPlasticity
end function plastic_nonlocal_init
2021-04-25 11:36:52 +05:30
module subroutine isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 13:09:17 +05:30
Mp
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine isotropic_LpAndItsTangent
2021-01-26 05:41:32 +05:30
2021-04-25 11:36:52 +05:30
pure module subroutine phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 13:09:17 +05:30
Mp
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine phenopowerlaw_LpAndItsTangent
2021-01-26 05:41:32 +05:30
2021-04-25 11:36:52 +05:30
pure module subroutine kinehardening_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 13:09:17 +05:30
Mp
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine kinehardening_LpAndItsTangent
2021-01-26 05:41:32 +05:30
module subroutine dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 13:09:17 +05:30
Mp
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine dislotwin_LpAndItsTangent
2021-01-26 05:41:32 +05:30
pure module subroutine dislotungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 13:09:17 +05:30
Mp
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine dislotungsten_LpAndItsTangent
2021-01-26 05:41:32 +05:30
module subroutine nonlocal_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
Lp
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2021-01-26 13:09:17 +05:30
dLp_dMp
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 05:41:32 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 12:03:04 +05:30
end subroutine nonlocal_LpAndItsTangent
2021-01-26 05:41:32 +05:30
2021-01-26 16:15:39 +05:30
module function isotropic_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
end function isotropic_dotState
2021-01-26 16:15:39 +05:30
module function phenopowerlaw_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
end function phenopowerlaw_dotState
2021-01-26 16:15:39 +05:30
module function plastic_kinehardening_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
end function plastic_kinehardening_dotState
2021-01-26 16:15:39 +05:30
2022-02-11 18:53:41 +05:30
module function dislotwin_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
2022-02-11 18:53:41 +05:30
dotState
end function dislotwin_dotState
2021-01-26 16:15:39 +05:30
module function dislotungsten_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
end function dislotungsten_dotState
2021-01-26 16:15:39 +05:30
2022-02-04 12:55:51 +05:30
module subroutine nonlocal_dotState(Mp,timestep,ph,en)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< MandelStress
real(pREAL), intent(in) :: &
2021-01-26 16:15:39 +05:30
timestep !< substepped crystallite time increment
integer, intent(in) :: &
ph, &
2022-02-04 12:55:51 +05:30
en
2021-01-26 16:47:00 +05:30
end subroutine nonlocal_dotState
2021-01-26 16:15:39 +05:30
2022-01-27 04:31:53 +05:30
module subroutine dislotwin_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 16:47:00 +05:30
end subroutine dislotwin_dependentState
2021-01-26 16:15:39 +05:30
2021-04-25 11:36:52 +05:30
module subroutine dislotungsten_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 16:47:00 +05:30
end subroutine dislotungsten_dependentState
2021-01-26 16:15:39 +05:30
2022-02-04 22:12:05 +05:30
module subroutine nonlocal_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2022-02-04 12:55:51 +05:30
en
2021-01-26 16:47:00 +05:30
end subroutine nonlocal_dependentState
2021-01-26 16:15:39 +05:30
2021-04-25 11:36:52 +05:30
module subroutine plastic_kinehardening_deltaState(Mp,ph,en)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 16:15:39 +05:30
end subroutine plastic_kinehardening_deltaState
2021-04-25 11:36:52 +05:30
module subroutine plastic_nonlocal_deltaState(Mp,ph,en)
real(pREAL), dimension(3,3), intent(in) :: &
2021-01-26 16:15:39 +05:30
Mp
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 16:15:39 +05:30
end subroutine plastic_nonlocal_deltaState
2021-01-26 05:41:32 +05:30
end interface
2021-01-26 13:09:17 +05:30
2021-01-26 05:41:32 +05:30
contains
2021-01-27 05:02:44 +05:30
module subroutine plastic_init
2021-01-27 15:14:03 +05:30
print'(/,1x,a)', '<<<+- phase:mechanical:plastic init -+>>>'
2021-01-27 15:14:03 +05:30
2021-12-11 14:24:46 +05:30
where(plastic_none_init()) phase_plasticity = PLASTIC_NONE_ID
where(plastic_isotropic_init()) phase_plasticity = PLASTIC_ISOTROPIC_ID
where(plastic_phenopowerlaw_init()) phase_plasticity = PLASTIC_PHENOPOWERLAW_ID
where(plastic_kinehardening_init()) phase_plasticity = PLASTIC_KINEHARDENING_ID
where(plastic_dislotwin_init()) phase_plasticity = PLASTIC_DISLOTWIN_ID
where(plastic_dislotungsten_init()) phase_plasticity = PLASTIC_DISLOTUNGSTEN_ID
where(plastic_nonlocal_init()) phase_plasticity = PLASTIC_NONLOCAL_ID
2021-01-27 05:02:44 +05:30
2021-12-11 14:24:46 +05:30
if (any(phase_plasticity == PLASTIC_undefined_ID)) call IO_error(201)
2021-01-27 05:02:44 +05:30
end subroutine plastic_init
2021-01-26 05:41:32 +05:30
!--------------------------------------------------------------------------------------------------
2023-02-11 20:00:12 +05:30
!> @brief constitutive equation for calculating the velocity gradient
2021-01-26 05:41:32 +05:30
! ToDo: Discuss whether it makes sense if crystallite handles the configuration conversion, i.e.
! Mp in, dLp_dMp out
!--------------------------------------------------------------------------------------------------
2021-01-26 16:11:19 +05:30
module subroutine plastic_LpAndItsTangents(Lp, dLp_dS, dLp_dFi, &
2021-04-13 00:51:15 +05:30
S, Fi, ph,en)
2021-01-26 05:41:32 +05:30
integer, intent(in) :: &
2021-04-13 00:51:15 +05:30
ph,en
real(pREAL), intent(in), dimension(3,3) :: &
2021-01-26 05:41:32 +05:30
S, & !< 2nd Piola-Kirchhoff stress
Fi !< intermediate deformation gradient
real(pREAL), intent(out), dimension(3,3) :: &
2021-01-26 05:41:32 +05:30
Lp !< plastic velocity gradient
real(pREAL), intent(out), dimension(3,3,3,3) :: &
2021-01-26 05:41:32 +05:30
dLp_dS, &
2021-04-13 00:51:15 +05:30
dLp_dFi !< derivative en Lp with respect to Fi
2021-01-26 05:41:32 +05:30
real(pREAL), dimension(3,3,3,3) :: &
2021-01-26 05:41:32 +05:30
dLp_dMp !< derivative of Lp with respect to Mandel stress
real(pREAL), dimension(3,3) :: &
2021-01-26 05:41:32 +05:30
Mp !< Mandel stress work conjugate with Lp
integer :: &
2021-02-14 21:59:23 +05:30
i, j
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
if (phase_plasticity(ph) == PLASTIC_NONE_ID) then
Lp = 0.0_pREAL
dLp_dFi = 0.0_pREAL
dLp_dS = 0.0_pREAL
else
2021-01-26 05:41:32 +05:30
Mp = matmul(matmul(transpose(Fi),Fi),S)
2021-02-14 21:59:23 +05:30
plasticType: select case (phase_plasticity(ph))
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_ISOTROPIC_ID) plasticType
call isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_PHENOPOWERLAW_ID) plasticType
call phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_KINEHARDENING_ID) plasticType
call kinehardening_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_NONLOCAL_ID) plasticType
call nonlocal_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTWIN_ID) plasticType
call dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTUNGSTEN_ID) plasticType
call dislotungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2021-01-26 05:41:32 +05:30
end select plasticType
2021-01-26 05:41:32 +05:30
do i=1,3; do j=1,3
dLp_dFi(i,j,1:3,1:3) = matmul(matmul(Fi,S),transpose(dLp_dMp(i,j,1:3,1:3))) + &
matmul(matmul(Fi,dLp_dMp(i,j,1:3,1:3)),S)
2022-02-01 13:53:03 +05:30
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
2021-01-26 05:41:32 +05:30
end if
2021-01-26 05:41:32 +05:30
2021-01-26 16:11:19 +05:30
end subroutine plastic_LpAndItsTangents
2021-01-26 13:09:17 +05:30
2021-01-26 16:15:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief contains the constitutive equation for calculating the rate of change of microstructure
!--------------------------------------------------------------------------------------------------
module function plastic_dotState(subdt,ph,en) result(dotState)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2021-04-13 00:51:15 +05:30
en
real(pREAL), intent(in) :: &
2021-01-26 16:15:39 +05:30
subdt !< timestep
real(pREAL), dimension(3,3) :: &
2021-01-26 16:15:39 +05:30
Mp
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
if (phase_plasticity(ph) /= PLASTIC_NONE_ID) then
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))
2021-01-26 16:15:39 +05:30
plasticType: select case (phase_plasticity(ph))
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_ISOTROPIC_ID) plasticType
dotState = isotropic_dotState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_PHENOPOWERLAW_ID) plasticType
dotState = phenopowerlaw_dotState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_KINEHARDENING_ID) plasticType
dotState = plastic_kinehardening_dotState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTWIN_ID) plasticType
2022-02-11 18:53:41 +05:30
dotState = dislotwin_dotState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTUNGSTEN_ID) plasticType
dotState = dislotungsten_dotState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_NONLOCAL_ID) plasticType
2022-02-04 12:55:51 +05:30
call nonlocal_dotState(Mp,subdt,ph,en)
dotState = plasticState(ph)%dotState(:,en)
end select plasticType
end if
2021-01-26 16:15:39 +05:30
2021-01-26 16:47:00 +05:30
end function plastic_dotState
2021-01-26 16:15:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calls microstructure function of the different plasticity constitutive models
!--------------------------------------------------------------------------------------------------
2022-02-04 22:12:05 +05:30
module subroutine plastic_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2021-01-26 16:15:39 +05:30
2021-05-22 14:56:19 +05:30
plasticType: select case (phase_plasticity(ph))
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTWIN_ID) plasticType
2022-01-27 04:31:53 +05:30
call dislotwin_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_DISLOTUNGSTEN_ID) plasticType
2021-04-25 11:36:52 +05:30
call dislotungsten_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_NONLOCAL_ID) plasticType
2022-02-04 12:55:51 +05:30
call nonlocal_dependentState(ph,en)
2021-01-26 16:15:39 +05:30
2021-02-09 03:51:53 +05:30
end select plasticType
2021-01-26 16:15:39 +05:30
2021-01-26 16:47:00 +05:30
end subroutine plastic_dependentState
2021-01-26 16:15:39 +05:30
!--------------------------------------------------------------------------------------------------
2023-02-11 20:00:12 +05:30
!> @brief for constitutive models that have an instantaneous change of state
2021-01-26 16:15:39 +05:30
!> will return false if delta state is not needed/supported by the constitutive model
!--------------------------------------------------------------------------------------------------
2021-04-13 00:51:15 +05:30
module function plastic_deltaState(ph, en) result(broken)
2021-01-26 16:15:39 +05:30
integer, intent(in) :: &
ph, &
2021-04-13 00:51:15 +05:30
en
logical :: broken
2021-01-26 16:15:39 +05:30
real(pREAL), dimension(3,3) :: &
2021-01-26 16:15:39 +05:30
Mp
integer :: &
mySize
2022-02-04 03:27:32 +05:30
broken = .false.
2021-01-26 16:15:39 +05:30
select case (phase_plasticity(ph))
2021-12-11 14:24:46 +05:30
case (PLASTIC_NONLOCAL_ID,PLASTIC_KINEHARDENING_ID)
2021-01-26 16:15:39 +05:30
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))
2021-01-26 16:15:39 +05:30
plasticType: select case (phase_plasticity(ph))
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_KINEHARDENING_ID) plasticType
call plastic_kinehardening_deltaState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
2021-12-11 14:24:46 +05:30
case (PLASTIC_NONLOCAL_ID) plasticType
call plastic_nonlocal_deltaState(Mp,ph,en)
2021-01-26 16:15:39 +05:30
end select plasticType
2021-01-26 16:15:39 +05:30
broken = any(IEEE_is_NaN(plasticState(ph)%deltaState(:,en)))
if (.not. broken) then
2021-01-26 16:15:39 +05:30
mySize = plasticState(ph)%sizeDeltaState
plasticState(ph)%deltaState2(1:mySize,en) = plasticState(ph)%deltaState2(1:mySize,en) &
+ plasticState(ph)%deltaState(1:mySize,en)
end if
end select
2021-01-26 16:15:39 +05:30
2021-01-26 16:47:00 +05:30
end function plastic_deltaState
2021-01-26 16:15:39 +05:30
2021-01-27 11:40:53 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief checks if a plastic module is active or not
!--------------------------------------------------------------------------------------------------
function plastic_active(plastic_label) result(active_plastic)
2021-01-27 11:40:53 +05:30
character(len=*), intent(in) :: plastic_label !< type of plasticity model
logical, dimension(:), allocatable :: active_plastic
type(tDict), pointer :: &
2021-01-27 11:40:53 +05:30
phases, &
phase, &
mech, &
pl
integer :: ph
phases => config_material%get_dict('phase')
2021-01-27 11:40:53 +05:30
allocate(active_plastic(phases%length), source = .false. )
do ph = 1, phases%length
phase => phases%get_dict(ph)
mech => phase%get_dict('mechanical')
pl => mech%get_dict('plastic',defaultVal = emptyDict)
2023-06-04 10:47:38 +05:30
active_plastic(ph) = pl%get_asStr('type',defaultVal='none') == plastic_label
2022-06-09 02:36:01 +05:30
end do
2021-01-27 11:40:53 +05:30
end function plastic_active
2021-01-26 05:41:32 +05:30
end submodule plastic