DAMASK_EICMD/src/phase_mechanical_plastic_ki...

477 lines
20 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) :: &
2020-09-22 18:56:34 +05:30
n = 1.0_pReal, & !< stress exponent for slip
dot_gamma_0 = 1.0_pReal !< reference shear strain rate for slip
real(pReal), allocatable, dimension(:) :: &
2020-09-22 18:56:34 +05:30
h_0_f, & !< initial hardening rate of forward stress for each slip
h_inf_f, & !< asymptotic hardening rate of forward stress for each slip
h_0_b, & !< initial hardening rate of back stress for each slip
h_inf_b, & !< asymptotic hardening rate of back stress for each slip
xi_inf_f, &
xi_inf_b
2021-07-21 18:33:28 +05:30
real(pReal), allocatable, dimension(:,:) :: &
h_sl_sl !< slip resistance from 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
2020-03-16 23:26:54 +05:30
logical :: &
nonSchmidActive = .false.
character(len=pStringLen), allocatable, dimension(:) :: &
output
character(len=:), allocatable, dimension(:) :: &
systems_sl
end type tParameters
type :: tIndexDotState
integer, dimension(2) :: &
xi, &
chi, &
gamma
end type tIndexDotState
type :: tKinehardeningState
2021-07-23 10:39:54 +05:30
real(pReal), pointer, dimension(:,:) :: &
xi, & !< resistance against plastic slip
chi, & !< back stress
chi_0, & !< back stress at last switch of stress sense
gamma, & !< accumulated (absolute) shear
gamma_0, & !< accumulated shear at last switch of stress sense
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, &
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
2020-03-16 23:26:54 +05:30
real(pReal), dimension(:), allocatable :: &
xi_0, & !< initial resistance against plastic flow
a !< non-Schmid coefficients
character(len=:), allocatable :: &
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 -+>>>'
print'(/,a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
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
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
#if defined (__GFORTRAN__)
2021-03-11 22:30:07 +05:30
prm%output = output_as1dString(pl)
2020-08-15 19:32:10 +05:30
#else
2021-03-11 22:30:07 +05:30
prm%output = pl%get_as1dString('output',defaultVal=emptyStringArray)
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
prm%systems_sl = lattice_labels_slip(N_sl,phase_lattice(ph))
2021-07-21 19:16:38 +05:30
prm%P = lattice_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph))
2021-07-21 19:16:38 +05:30
if (phase_lattice(ph) == 'cI') then
2021-03-11 22:30:07 +05:30
a = pl%get_as1dFloat('a_nonSchmid',defaultVal = emptyRealArray)
2022-12-07 22:59:03 +05:30
if (size(a) > 0) prm%nonSchmidActive = .true.
2021-07-21 19:59:57 +05:30
prm%P_nS_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
prm%P_nS_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
else
2021-07-21 19:59:57 +05:30
prm%P_nS_pos = prm%P
prm%P_nS_neg = prm%P
end if
2021-07-21 19:16:38 +05:30
prm%h_sl_sl = lattice_interaction_SlipBySlip(N_sl,pl%get_as1dFloat('h_sl-sl'), &
phase_lattice(ph))
2021-03-11 22:30:07 +05:30
xi_0 = pl%get_as1dFloat('xi_0', requiredSize=size(N_sl))
prm%xi_inf_f = pl%get_as1dFloat('xi_inf_f', requiredSize=size(N_sl))
prm%xi_inf_b = pl%get_as1dFloat('xi_inf_b', requiredSize=size(N_sl))
prm%h_0_f = pl%get_as1dFloat('h_0_f', requiredSize=size(N_sl))
prm%h_inf_f = pl%get_as1dFloat('h_inf_f', requiredSize=size(N_sl))
prm%h_0_b = pl%get_as1dFloat('h_0_b', requiredSize=size(N_sl))
prm%h_inf_b = pl%get_as1dFloat('h_inf_b', requiredSize=size(N_sl))
2020-09-22 18:56:34 +05:30
prm%dot_gamma_0 = pl%get_asFloat('dot_gamma_0')
prm%n = pl%get_asFloat('n')
! expand: family => system
2020-09-22 18:56:34 +05:30
xi_0 = math_expand(xi_0, N_sl)
prm%xi_inf_f = math_expand(prm%xi_inf_f, N_sl)
prm%xi_inf_b = math_expand(prm%xi_inf_b, N_sl)
prm%h_0_f = math_expand(prm%h_0_f, N_sl)
prm%h_inf_f = math_expand(prm%h_inf_f, N_sl)
prm%h_0_b = math_expand(prm%h_0_b, N_sl)
prm%h_inf_b = math_expand(prm%h_inf_b, N_sl)
!--------------------------------------------------------------------------------------------------
! sanity checks
2020-09-22 18:56:34 +05:30
if ( prm%dot_gamma_0 <= 0.0_pReal) extmsg = trim(extmsg)//' dot_gamma_0'
if ( 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_f <= 0.0_pReal)) extmsg = trim(extmsg)//' xi_inf_f'
if (any(prm%xi_inf_b <= 0.0_pReal)) extmsg = trim(extmsg)//' xi_inf_b'
2020-03-16 23:26:54 +05:30
else slipActive
xi_0 = emptyRealArray
2020-09-22 18:56:34 +05:30
allocate(prm%xi_inf_f,prm%xi_inf_b,prm%h_0_f,prm%h_inf_f,prm%h_0_b,prm%h_inf_b,source=emptyRealArray)
2021-07-21 18:33:28 +05:30
allocate(prm%h_sl_sl(0,0))
end if slipActive
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2021-04-06 15:08:44 +05:30
Nmembers = count(material_phaseID == ph)
2021-07-23 10:39:54 +05:30
sizeDotState = size(['xi ','chi ', 'gamma']) * prm%sum_N_sl
sizeDeltaState = size(['sgn_gamma', 'chi_0 ', 'gamma_0 ']) * prm%sum_N_sl
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)
2021-02-16 20:36:09 +05:30
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_xi',defaultVal=1.0_pReal)
2022-12-07 22:59:03 +05:30
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,:)
2021-02-16 20:36:09 +05:30
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('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,:)
2021-02-16 20:36:09 +05:30
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_gamma',defaultVal=1.0e-6_pReal)
2022-12-07 22:59:03 +05:30
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
2021-07-23 10:39:54 +05:30
stt%chi_0 => plasticState(ph)%state (startIndex :endIndex ,:)
dlt%chi_0 => plasticState(ph)%deltaState(startIndex-o:endIndex-o,:)
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%gamma_0 => plasticState(ph)%state (startIndex :endIndex ,:)
dlt%gamma_0 => 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.
!--------------------------------------------------------------------------------------------------
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) :: &
Lp !< plastic velocity gradient
real(pReal), dimension(3,3,3,3), intent(out) :: &
dLp_dMp !< derivative of Lp with respect to the Mandel stress
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
integer :: &
i,k,l,m,n
2021-02-14 05:20:42 +05:30
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_pos,dot_gamma_neg, &
ddot_gamma_dtau_pos,ddot_gamma_dtau_neg
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
2021-02-14 05:20:42 +05:30
associate(prm => param(ph))
2021-07-21 20:04:33 +05:30
call kinetics(Mp,ph,en,dot_gamma_pos,dot_gamma_neg,ddot_gamma_dtau_pos,ddot_gamma_dtau_neg)
2020-03-16 01:46:28 +05:30
do i = 1, prm%sum_N_sl
2021-07-21 20:04:33 +05:30
Lp = Lp + (dot_gamma_pos(i)+dot_gamma_neg(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) &
2021-07-21 20:04:33 +05:30
+ ddot_gamma_dtau_pos(i) * prm%P(k,l,i) * prm%P_nS_pos(m,n,i) &
+ ddot_gamma_dtau_neg(i) * prm%P(k,l,i) * prm%P_nS_neg(m,n,i)
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
2021-02-14 05:20:42 +05:30
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_pos,dot_gamma_neg
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)))
2021-07-23 03:39:51 +05:30
call kinetics(Mp,ph,en,dot_gamma_pos,dot_gamma_neg)
dot_gamma = abs(dot_gamma_pos+dot_gamma_neg)
2021-07-23 10:39:54 +05:30
sumGamma = sum(stt%gamma(:,en))
dot_xi = matmul(prm%h_sl_sl,dot_gamma) &
2021-07-23 03:39:51 +05:30
* ( prm%h_inf_f &
+ (prm%h_0_f - prm%h_inf_f + prm%h_0_f*prm%h_inf_f*sumGamma/prm%xi_inf_f) &
* exp(-sumGamma*prm%h_0_f/prm%xi_inf_f) &
)
dot_chi = stt%sgn_gamma(:,en)*dot_gamma &
* ( prm%h_inf_b &
+ (prm%h_0_b - prm%h_inf_b &
2021-07-23 10:39:54 +05:30
+ prm%h_0_b*prm%h_inf_b/(prm%xi_inf_b+stt%chi_0(:,en))*(stt%gamma(:,en)-stt%gamma_0(:,en))&
) *exp(-(stt%gamma(:,en)-stt%gamma_0(:,en)) *prm%h_0_b/(prm%xi_inf_b+stt%chi_0(:,en))) &
2021-07-23 03:39:51 +05:30
)
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
2021-02-14 05:20:42 +05:30
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_pos,dot_gamma_neg, &
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))
2021-07-23 10:39:54 +05:30
call kinetics(Mp,ph,en,dot_gamma_pos,dot_gamma_neg)
sgn_gamma = merge(state(ph)%sgn_gamma(:,en), &
sign(1.0_pReal,dot_gamma_pos+dot_gamma_neg), &
dEq0(dot_gamma_pos+dot_gamma_neg,1e-10_pReal))
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_0 (:,en) = abs(stt%chi(:,en)) - stt%chi_0(:,en)
dlt%gamma_0(:,en) = stt%gamma(:,en) - stt%gamma_0(:,en)
else where
dlt%sgn_gamma (:,en) = 0.0_pReal
dlt%chi_0 (:,en) = 0.0_pReal
dlt%gamma_0(:,en) = 0.0_pReal
end where
end associate
end subroutine plastic_kinehardening_deltaState
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2021-02-14 05:20:42 +05:30
module subroutine plastic_kinehardening_results(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')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%xi,group,trim(prm%output(ou)), &
'resistance against plastic slip','Pa',prm%systems_sl)
case ('chi')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%chi,group,trim(prm%output(ou)), &
'back stress','Pa',prm%systems_sl)
case ('sgn(gamma)')
2021-08-08 01:32:44 +05:30
call results_writeDataset(int(stt%sgn_gamma),group,trim(prm%output(ou)), &
'sense of shear','1',prm%systems_sl)
case ('chi_0')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%chi_0,group,trim(prm%output(ou)), &
'back stress at last switch of stress sense','Pa',prm%systems_sl)
case ('gamma_0')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%gamma_0,group,trim(prm%output(ou)), &
'plastic shear at last switch of stress sense','1',prm%systems_sl)
case ('gamma')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%gamma,group,trim(prm%output(ou)), &
'plastic shear','1',prm%systems_sl)
end select
end do
end associate
end subroutine plastic_kinehardening_results
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.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end.
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
pure subroutine kinetics(Mp,ph,en, &
2021-07-21 20:04:33 +05:30
dot_gamma_pos,dot_gamma_neg,ddot_gamma_dtau_pos,ddot_gamma_dtau_neg)
2018-11-26 11:40:43 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2021-07-23 10:39:54 +05:30
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
2021-02-14 05:20:42 +05:30
real(pReal), intent(out), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_pos, &
dot_gamma_neg
2021-07-23 10:39:54 +05:30
real(pReal), intent(out), dimension(param(ph)%sum_N_sl), optional :: &
2021-07-21 20:04:33 +05:30
ddot_gamma_dtau_pos, &
ddot_gamma_dtau_neg
2021-02-14 05:20:42 +05:30
real(pReal), dimension(param(ph)%sum_N_sl) :: &
tau_pos, &
tau_neg
integer :: i
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-23 03:39:51 +05:30
do i = 1, prm%sum_N_sl
2021-07-23 10:39:54 +05:30
tau_pos(i) = math_tensordot(Mp,prm%P_nS_pos(1:3,1:3,i)) - stt%chi(i,en)
tau_neg(i) = merge(math_tensordot(Mp,prm%P_nS_neg(1:3,1:3,i)) - stt%chi(i,en), &
2021-07-23 03:39:51 +05:30
0.0_pReal, prm%nonSchmidActive)
end do
2021-07-23 03:39:51 +05:30
where(dNeq0(tau_pos))
dot_gamma_pos = prm%dot_gamma_0 * merge(0.5_pReal,1.0_pReal, prm%nonSchmidActive) & ! 1/2 if non-Schmid active
2021-07-23 10:39:54 +05:30
* sign(abs(tau_pos/stt%xi(:,en))**prm%n, tau_pos)
else where
2021-07-23 03:39:51 +05:30
dot_gamma_pos = 0.0_pReal
end where
2021-07-23 03:39:51 +05:30
where(dNeq0(tau_neg))
dot_gamma_neg = prm%dot_gamma_0 * 0.5_pReal & ! only used if non-Schmid active, always 1/2
2021-07-23 10:39:54 +05:30
* sign(abs(tau_neg/stt%xi(:,en))**prm%n, tau_neg)
else where
2021-07-23 03:39:51 +05:30
dot_gamma_neg = 0.0_pReal
end where
2021-07-23 03:39:51 +05:30
if (present(ddot_gamma_dtau_pos)) then
where(dNeq0(dot_gamma_pos))
ddot_gamma_dtau_pos = dot_gamma_pos*prm%n/tau_pos
else where
ddot_gamma_dtau_pos = 0.0_pReal
end where
end if
2021-07-23 03:39:51 +05:30
if (present(ddot_gamma_dtau_neg)) then
where(dNeq0(dot_gamma_neg))
ddot_gamma_dtau_neg = dot_gamma_neg*prm%n/tau_neg
else where
ddot_gamma_dtau_neg = 0.0_pReal
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