DAMASK_EICMD/src/phase_mechanical_plastic_ki...

475 lines
21 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Philip Eisenlohr, Michigan State University
!> @author Zhuowen Zhao, Michigan State University
2018-11-25 15:44:09 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Phenomenological crystal plasticity using a power-law formulation for the shear rates
!! and a Voce-type kinematic hardening rule.
!--------------------------------------------------------------------------------------------------
submodule(phase:plastic) kinehardening
type :: tParameters
real(pREAL), allocatable, dimension(:) :: &
dot_gamma_0, & !< reference shear strain rate for slip
n, & !< stress exponent for slip
h_0_xi, & !< initial hardening rate of forest stress per slip family
!! θ_0,for
h_0_chi, & !< initial hardening rate of back stress per slip family
!! θ_0,bs
h_inf_xi, & !< asymptotic hardening rate of forest stress per slip family
!! θ_1,for
h_inf_chi, & !< asymptotic hardening rate of back stress per slip family
!! θ_1,bs
xi_inf, & !< back-extrapolated forest stress from terminal linear hardening
chi_inf !< back-extrapolated back stress from terminal linear hardening
real(pREAL), allocatable, dimension(:,:) :: &
h_sl_sl !< slip resistance change per slip activity
real(pREAL), allocatable, dimension(:,:,:) :: &
2020-03-16 01:46:28 +05:30
P, &
2021-07-21 19:59:57 +05:30
P_nS_pos, &
P_nS_neg
integer :: &
2021-02-16 20:36:09 +05:30
sum_N_sl
2023-06-04 10:47:38 +05:30
character(len=pSTRLEN), allocatable, dimension(:) :: &
output
2023-06-04 10:47:38 +05:30
character(len=:), allocatable, dimension(:) :: &
systems_sl
end type tParameters
type :: tIndexDotState
integer, dimension(2) :: &
xi, &
chi, &
gamma
end type tIndexDotState
type :: tKinehardeningState
real(pREAL), pointer, dimension(:,:) :: &
xi, & !< forest stress
!! τ_for
2021-07-23 10:39:54 +05:30
chi, & !< back stress
!! τ_bs
chi_flip, & !< back stress at last reversal of stress sense
!! χ_0
2021-07-23 10:39:54 +05:30
gamma, & !< accumulated (absolute) shear
gamma_flip, & !< accumulated shear at last reversal of stress sense
!! γ_0
2021-07-23 10:39:54 +05:30
sgn_gamma !< sense of acting shear stress (-1 or +1)
end type tKinehardeningState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:) :: param
type(tIndexDotState), allocatable, dimension(:) :: indexDotState
type(tKinehardeningState), allocatable, dimension(:) :: state, deltaState
contains
!--------------------------------------------------------------------------------------------------
!> @brief Perform module initialization.
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
module function plastic_kinehardening_init() result(myPlasticity)
2020-08-15 19:32:10 +05:30
logical, dimension(:), allocatable :: myPlasticity
integer :: &
2021-02-16 20:36:09 +05:30
ph, o, &
N_fam, &
2021-03-05 01:46:36 +05:30
Nmembers, &
sizeState, sizeDeltaState, sizeDotState, &
2020-07-02 00:52:05 +05:30
startIndex, endIndex
2020-03-16 23:26:54 +05:30
integer, dimension(:), allocatable :: &
2020-03-16 03:37:23 +05:30
N_sl
real(pREAL), dimension(:), allocatable :: &
2023-09-15 20:43:58 +05:30
xi_0 !< initial forest stress
!! τ_for,0
2023-09-15 20:43:58 +05:30
real(pREAL), dimension(:,:), allocatable :: &
a_nS !< non-Schmid coefficients
character(len=:), allocatable :: &
refs, &
extmsg
type(tDict), pointer :: &
2020-08-15 19:32:10 +05:30
phases, &
phase, &
2020-11-03 03:16:46 +05:30
mech, &
2020-08-15 19:32:10 +05:30
pl
2020-08-15 19:32:10 +05:30
myPlasticity = plastic_active('kinehardening')
2022-12-07 22:59:03 +05:30
if (count(myPlasticity) == 0) return
2020-07-02 00:52:05 +05:30
print'(/,1x,a)', '<<<+- phase:mechanical:plastic:kinehardening init -+>>>'
2021-02-13 17:37:35 +05:30
print'(/,1x,a)', 'J.A. Wollmershauser et al., International Journal of Fatigue 36:181193, 2012'
print'( 1x,a)', 'https://doi.org/10.1016/j.ijfatigue.2011.07.008'
2021-02-13 17:37:35 +05:30
2023-08-29 07:11:48 +05:30
print'(/,1x,a,1x,i0)', '# phases:',count(myPlasticity); flush(IO_STDOUT)
2023-01-18 23:20:01 +05:30
phases => config_material%get_dict('phase')
2021-02-14 05:20:42 +05:30
allocate(param(phases%length))
allocate(indexDotState(phases%length))
2021-02-14 05:20:42 +05:30
allocate(state(phases%length))
allocate(deltaState(phases%length))
extmsg = ''
2021-02-14 05:20:42 +05:30
2021-02-16 20:36:09 +05:30
do ph = 1, phases%length
if (.not. myPlasticity(ph)) cycle
2021-02-16 20:36:09 +05:30
associate(prm => param(ph), &
stt => state(ph), dlt => deltaState(ph), &
idx_dot => indexDotState(ph))
2021-02-16 20:36:09 +05:30
phase => phases%get_dict(ph)
mech => phase%get_dict('mechanical')
pl => mech%get_dict('plastic')
2020-08-15 19:32:10 +05:30
print'(/,1x,a,1x,i0,a)', 'phase',ph,': '//phases%key(ph)
refs = config_listReferences(pl,indent=3)
if (len(refs) > 0) print'(/,1x,a)', refs
2020-08-15 19:32:10 +05:30
#if defined (__GFORTRAN__)
2023-06-04 10:47:38 +05:30
prm%output = output_as1dStr(pl)
2020-08-15 19:32:10 +05:30
#else
2023-06-04 10:47:38 +05:30
prm%output = pl%get_as1dStr('output',defaultVal=emptyStrArray)
2020-08-15 19:32:10 +05:30
#endif
2020-03-15 00:33:57 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
2021-03-11 22:30:07 +05:30
N_sl = pl%get_as1dInt('N_sl',defaultVal=emptyIntArray)
2020-03-16 23:26:54 +05:30
prm%sum_N_sl = sum(abs(N_sl))
2020-03-16 01:46:28 +05:30
slipActive: if (prm%sum_N_sl > 0) then
N_fam = size(N_sl)
2023-07-15 23:58:57 +05:30
prm%systems_sl = crystal_labels_slip(N_sl,phase_lattice(ph))
prm%P = crystal_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph))
2021-07-21 19:16:38 +05:30
if (phase_lattice(ph) == 'cI') then
2023-09-15 20:43:58 +05:30
allocate(a_nS(3,size(pl%get_as1dReal('a_nonSchmid_110',defaultVal=emptyRealArray))),source=0.0_pREAL)
a_nS(1,:) = pl%get_as1dReal('a_nonSchmid_110',defaultVal=emptyRealArray)
prm%P_nS_pos = crystal_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph),nonSchmidCoefficients=a_nS,sense=+1)
prm%P_nS_neg = crystal_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph),nonSchmidCoefficients=a_nS,sense=-1)
else
2023-09-15 20:43:58 +05:30
prm%P_nS_pos = +prm%P
prm%P_nS_neg = -prm%P
end if
2023-07-15 23:58:57 +05:30
prm%h_sl_sl = crystal_interaction_SlipBySlip(N_sl,pl%get_as1dReal('h_sl-sl'), &
phase_lattice(ph))
xi_0 = math_expand(pl%get_as1dReal('xi_0', requiredSize=N_fam),N_sl)
prm%dot_gamma_0 = math_expand(pl%get_as1dReal('dot_gamma_0', requiredSize=N_fam),N_sl)
prm%n = math_expand(pl%get_as1dReal('n', requiredSize=N_fam),N_sl)
prm%xi_inf = math_expand(pl%get_as1dReal('xi_inf', requiredSize=N_fam),N_sl)
prm%chi_inf = math_expand(pl%get_as1dReal('chi_inf', requiredSize=N_fam),N_sl)
prm%h_0_xi = math_expand(pl%get_as1dReal('h_0_xi', requiredSize=N_fam),N_sl)
prm%h_0_chi = math_expand(pl%get_as1dReal('h_0_chi', requiredSize=N_fam),N_sl)
prm%h_inf_xi = math_expand(pl%get_as1dReal('h_inf_xi', requiredSize=N_fam),N_sl)
prm%h_inf_chi = math_expand(pl%get_as1dReal('h_inf_chi', requiredSize=N_fam),N_sl)
!--------------------------------------------------------------------------------------------------
! sanity checks
if (any(prm%dot_gamma_0 <= 0.0_pREAL)) extmsg = trim(extmsg)//' dot_gamma_0'
if (any(prm%n <= 0.0_pREAL)) extmsg = trim(extmsg)//' n'
if (any(xi_0 <= 0.0_pREAL)) extmsg = trim(extmsg)//' xi_0'
if (any(prm%xi_inf <= 0.0_pREAL)) extmsg = trim(extmsg)//' xi_inf'
if (any(prm%chi_inf <= 0.0_pREAL)) extmsg = trim(extmsg)//' chi_inf'
2020-03-16 23:26:54 +05:30
else slipActive
xi_0 = emptyRealArray
allocate(prm%dot_gamma_0, &
prm%n, &
prm%xi_inf, &
prm%chi_inf, &
prm%h_0_xi, &
prm%h_0_chi, &
prm%h_inf_xi, &
prm%h_inf_chi, &
source=emptyRealArray)
2021-07-21 18:33:28 +05:30
allocate(prm%h_sl_sl(0,0))
end if slipActive
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2023-01-23 13:01:59 +05:30
Nmembers = count(material_ID_phase == ph)
sizeDotState = prm%sum_N_sl * size(['xi ',&
'chi ',&
'gamma'])
sizeDeltaState = prm%sum_N_sl * size(['sgn_gamma ',&
'chi_flip ',&
'gamma_flip'])
sizeState = sizeDotState + sizeDeltaState
2021-03-05 01:46:36 +05:30
call phase_allocateState(plasticState(ph),Nmembers,sizeState,sizeDotState,sizeDeltaState)
deallocate(plasticState(ph)%dotState) ! ToDo: remove dotState completely
!--------------------------------------------------------------------------------------------------
2020-03-16 03:37:23 +05:30
! state aliases and initialization
startIndex = 1
2020-03-16 01:46:28 +05:30
endIndex = prm%sum_N_sl
idx_dot%xi = [startIndex,endIndex]
stt%xi => plasticState(ph)%state(startIndex:endIndex,:)
2021-07-23 10:39:54 +05:30
stt%xi = spread(xi_0, 2, Nmembers)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_xi',defaultVal=1.0_pREAL)
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_xi'
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
idx_dot%chi = [startIndex,endIndex]
stt%chi => plasticState(ph)%state(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_xi',defaultVal=1.0_pREAL)
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
idx_dot%gamma = [startIndex,endIndex]
stt%gamma => plasticState(ph)%state(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_gamma',defaultVal=1.0e-6_pREAL)
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_gamma'
2021-02-16 20:36:09 +05:30
o = plasticState(ph)%offsetDeltaState
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
2021-07-23 10:39:54 +05:30
stt%sgn_gamma => plasticState(ph)%state (startIndex :endIndex ,:)
dlt%sgn_gamma => plasticState(ph)%deltaState(startIndex-o:endIndex-o,:)
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
stt%chi_flip => plasticState(ph)%state (startIndex :endIndex ,:)
dlt%chi_flip => plasticState(ph)%deltaState(startIndex-o:endIndex-o,:)
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
stt%gamma_flip => plasticState(ph)%state (startIndex :endIndex ,:)
dlt%gamma_flip => plasticState(ph)%deltaState(startIndex-o:endIndex-o,:)
end associate
2020-03-15 00:33:57 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg))
2020-03-15 00:33:57 +05:30
end do
2020-03-15 00:33:57 +05:30
2020-08-15 19:32:10 +05:30
end function plastic_kinehardening_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
2023-09-15 20:43:58 +05:30
pure module subroutine kinehardening_LpAndItsTangent(Lp,dLp_dMp, Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
Lp !< plastic velocity gradient
real(pREAL), dimension(3,3,3,3), intent(out) :: &
dLp_dMp !< derivative of Lp with respect to the Mandel stress
2023-09-15 20:43:58 +05:30
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2023-09-15 20:43:58 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
integer :: &
i,k,l,m,n
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2023-09-15 20:43:58 +05:30
dot_gamma, ddot_gamma_dtau
real(pREAL), dimension(3,3,param(ph)%sum_N_sl) :: &
P_nS
2023-09-22 10:40:08 +05:30
Lp = 0.0_pREAL
dLp_dMp = 0.0_pREAL
2021-02-14 05:20:42 +05:30
associate(prm => param(ph))
2023-09-22 10:40:08 +05:30
call kinetics(Mp,ph,en, dot_gamma,ddot_gamma_dtau)
P_nS = merge(prm%P_nS_pos,prm%P_nS_neg, spread(spread(dot_gamma,1,3),2,3)>0.0_pREAL) ! faster than 'merge' in loop
2023-09-22 10:40:08 +05:30
do i = 1, prm%sum_N_sl
Lp = Lp + dot_gamma(i)*prm%P(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
+ ddot_gamma_dtau(i) * prm%P(k,l,i) * P_nS(m,n,i)
2023-09-22 10:40:08 +05:30
end do
end associate
2021-01-26 12:03:04 +05:30
end subroutine kinehardening_LpAndItsTangent
2018-12-22 03:11:39 +05:30
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate the rate of change of microstructure.
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
module function plastic_kinehardening_dotState(Mp,ph,en) result(dotState)
2019-01-07 12:06:11 +05:30
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
real(pREAL) :: &
sumGamma
associate(prm => param(ph), stt => state(ph), &
dot_xi => dotState(IndexDotState(ph)%xi(1):IndexDotState(ph)%xi(2)),&
dot_chi => dotState(IndexDotState(ph)%chi(1):IndexDotState(ph)%chi(2)),&
dot_gamma => dotState(IndexDotState(ph)%gamma(1):IndexDotState(ph)%gamma(2)))
2023-09-15 20:43:58 +05:30
call kinetics(Mp,ph,en, dot_gamma)
2021-07-23 10:39:54 +05:30
sumGamma = sum(stt%gamma(:,en))
2023-09-15 20:43:58 +05:30
dot_gamma = abs(dot_gamma)
dot_xi = matmul(prm%h_sl_sl,dot_gamma) &
* ( prm%h_inf_xi &
+ ( prm%h_0_xi &
- prm%h_inf_xi * (1_pREAL -sumGamma*prm%h_0_xi/prm%xi_inf) ) &
* exp(-sumGamma*prm%h_0_xi/prm%xi_inf) &
)
dot_chi = stt%sgn_gamma(:,en)*dot_gamma &
* ( prm%h_inf_chi &
+ ( prm%h_0_chi &
- prm%h_inf_chi*(1_pREAL -(stt%gamma(:,en)-stt%gamma_flip(:,en))*prm%h_0_chi/(prm%chi_inf+stt%chi_flip(:,en))) ) &
* exp(-(stt%gamma(:,en)-stt%gamma_flip(:,en))*prm%h_0_chi/(prm%chi_inf+stt%chi_flip(:,en))) &
)
end associate
2019-01-07 12:06:11 +05:30
end function plastic_kinehardening_dotState
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate (instantaneous) incremental change of microstructure.
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine plastic_kinehardening_deltaState(Mp,ph,en)
2018-12-30 20:09:48 +05:30
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2023-09-15 20:43:58 +05:30
dot_gamma, &
2021-07-23 10:39:54 +05:30
sgn_gamma
2021-07-23 10:39:54 +05:30
associate(prm => param(ph), stt => state(ph), dlt => deltaState(ph))
2023-09-15 20:43:58 +05:30
call kinetics(Mp,ph,en, dot_gamma)
2021-07-23 10:39:54 +05:30
sgn_gamma = merge(state(ph)%sgn_gamma(:,en), &
2023-09-15 20:43:58 +05:30
sign(1.0_pREAL,dot_gamma), &
dEq0(dot_gamma,1e-10_pREAL))
2021-07-23 10:39:54 +05:30
where(dNeq(sgn_gamma,stt%sgn_gamma(:,en),0.1_pREAL)) ! ToDo sgn_gamma*stt%sgn_gamma(:,en)<0
dlt%sgn_gamma (:,en) = sgn_gamma - stt%sgn_gamma (:,en)
dlt%chi_flip (:,en) = abs(stt%chi (:,en)) - stt%chi_flip (:,en)
dlt%gamma_flip(:,en) = stt%gamma(:,en) - stt%gamma_flip(:,en)
2021-07-23 10:39:54 +05:30
else where
dlt%sgn_gamma (:,en) = 0.0_pREAL
dlt%chi_flip (:,en) = 0.0_pREAL
dlt%gamma_flip(:,en) = 0.0_pREAL
2021-07-23 10:39:54 +05:30
end where
end associate
end subroutine plastic_kinehardening_deltaState
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
module subroutine plastic_kinehardening_result(ph,group)
2021-02-14 05:20:42 +05:30
integer, intent(in) :: ph
character(len=*), intent(in) :: group
integer :: ou
2021-02-14 05:20:42 +05:30
associate(prm => param(ph), stt => state(ph))
do ou = 1,size(prm%output)
select case(trim(prm%output(ou)))
case ('xi')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%xi,group,trim(prm%output(ou)), &
'forest stress','Pa',prm%systems_sl)
case ('chi')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%chi,group,trim(prm%output(ou)), &
'back stress','Pa',prm%systems_sl)
case ('sgn(gamma)')
2023-01-19 22:07:45 +05:30
call result_writeDataset(int(stt%sgn_gamma),group,trim(prm%output(ou)), &
'sense of shear','1',prm%systems_sl)
case ('chi_flip')
call result_writeDataset(stt%chi_flip,group,trim(prm%output(ou)), &
'back stress at last reversal of stress sense','Pa',prm%systems_sl)
case ('gamma_flip')
call result_writeDataset(stt%gamma_flip,group,trim(prm%output(ou)), &
'plastic shear at last reversal of stress sense','1',prm%systems_sl)
case ('gamma')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%gamma,group,trim(prm%output(ou)), &
'plastic shear','1',prm%systems_sl)
end select
end do
end associate
2023-01-19 22:07:45 +05:30
end subroutine plastic_kinehardening_result
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate shear rates on slip systems and their derivatives with respect to resolved
! stress.
!> @details: Derivatives are calculated only optionally.
2023-02-11 20:00:12 +05:30
! NOTE: Contrary to common convention, here the result (i.e. intent(out)) variables have to be put
! at the end since some of them are optional.
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
pure subroutine kinetics(Mp,ph,en, &
2023-09-15 20:43:58 +05:30
dot_gamma,ddot_gamma_dtau)
2018-11-26 11:40:43 +05:30
2023-09-15 20:43:58 +05:30
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2023-09-15 20:43:58 +05:30
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
2023-09-15 20:43:58 +05:30
real(pREAL), dimension(param(ph)%sum_N_sl), intent(out) :: &
dot_gamma
real(pREAL), dimension(param(ph)%sum_N_sl), optional, intent(out) :: &
ddot_gamma_dtau
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
tau_pos, &
tau_neg
integer :: i
2023-09-22 10:40:08 +05:30
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph))
2023-09-15 20:43:58 +05:30
tau_pos = [(math_tensordot(Mp,prm%P_nS_pos(1:3,1:3,i)) - stt%chi(i,en),i=1,prm%sum_N_sl)]
tau_neg = [(math_tensordot(Mp,prm%P_nS_neg(1:3,1:3,i)) + stt%chi(i,en),i=1,prm%sum_N_sl)]
2021-07-23 03:39:51 +05:30
2023-09-15 20:43:58 +05:30
dot_gamma = merge(+1.0_pREAL,-1.0_pREAL, tau_pos>tau_neg) &
* prm%dot_gamma_0 &
2023-09-22 02:02:38 +05:30
* (max(tau_pos,tau_neg)/stt%xi(:,en))**prm%n
2021-07-23 03:39:51 +05:30
2023-09-15 20:43:58 +05:30
if (present(ddot_gamma_dtau)) then
where(dNeq0(dot_gamma))
2023-09-22 02:02:38 +05:30
ddot_gamma_dtau = dot_gamma*prm%n/max(tau_pos,tau_neg)
2021-07-23 03:39:51 +05:30
else where
2023-09-15 20:43:58 +05:30
ddot_gamma_dtau = 0.0_pREAL
2021-07-23 03:39:51 +05:30
end where
end if
2021-07-23 03:39:51 +05:30
end associate
2018-11-26 11:40:43 +05:30
end subroutine kinetics
2021-01-26 05:41:32 +05:30
end submodule kinehardening