DAMASK_EICMD/src/phase_mechanical_plastic_di...

531 lines
24 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @author David Cereceda, Lawrence Livermore National Laboratory
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief crystal plasticity model for bcc metals, especially Tungsten
!--------------------------------------------------------------------------------------------------
submodule(phase:plastic) dislotungsten
2019-12-03 02:21:25 +05:30
type :: tParameters
real(pREAL) :: &
D = 1.0_pREAL, & !< grain size
D_0 = 1.0_pREAL, & !< prefactor for self-diffusion coefficient
Q_cl = 1.0_pREAL !< activation energy for dislocation climb
real(pREAL), allocatable, dimension(:) :: &
2021-04-26 13:56:50 +05:30
b_sl, & !< magnitude of Burgers vector [m]
2021-07-24 01:34:39 +05:30
d_caron, & !< distance of spontaneous annhihilation
2019-03-22 15:25:17 +05:30
i_sl, & !< Adj. parameter for distance between 2 forest dislocations
2020-09-22 18:05:05 +05:30
f_at, & !< factor to calculate atomic volume
2020-09-23 04:25:19 +05:30
tau_Peierls, & !< Peierls stress
2019-03-22 15:25:17 +05:30
!* mobility law parameters
2020-09-22 18:05:05 +05:30
Q_s, & !< activation energy for glide [J]
2019-03-22 15:25:17 +05:30
p, & !< p-exponent in glide velocity
q, & !< q-exponent in glide velocity
B, & !< friction coefficient
2020-09-22 18:05:05 +05:30
h, & !< height of the kink pair
2019-03-22 15:25:17 +05:30
w, & !< width of the kink pair
omega !< attempt frequency for kink pair nucleation
real(pREAL), allocatable, dimension(:,:) :: &
2019-03-22 15:25:17 +05:30
h_sl_sl, & !< slip resistance from slip activity
2020-03-16 16:23:36 +05:30
forestProjection
real(pREAL), allocatable, dimension(:,:,:) :: &
2020-03-16 18:03:14 +05:30
P_sl, &
2021-07-21 19:59:57 +05:30
P_nS_pos, &
P_nS_neg
2019-03-22 15:25:17 +05:30
integer :: &
sum_N_sl !< total number of active slip system
2022-11-30 18:03:04 +05:30
character(len=:), allocatable :: &
isotropic_bound
2023-06-04 10:47:38 +05:30
character(len=pSTRLEN), allocatable, dimension(:) :: &
output
2019-03-22 15:25:17 +05:30
logical :: &
dipoleFormation !< flag indicating consideration of dipole formation
character(len=:), allocatable, dimension(:) :: &
systems_sl
end type tParameters !< container type for internal constitutive parameters
type :: tIndexDotState
integer, dimension(2) :: &
rho_mob, &
rho_dip, &
gamma_sl
end type tIndexDotState
type :: tDislotungstenState
real(pREAL), dimension(:,:), pointer :: &
2019-03-22 15:25:17 +05:30
rho_mob, &
rho_dip, &
gamma_sl
end type tDislotungstenState
2021-07-27 02:26:40 +05:30
type :: tDislotungstenDependentState
real(pREAL), dimension(:,:), allocatable :: &
2019-03-22 15:25:17 +05:30
Lambda_sl, &
2021-07-24 01:26:34 +05:30
tau_pass
2021-07-27 02:26:40 +05:30
end type tDislotungstenDependentState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:) :: param
type(tIndexDotState), allocatable, dimension(:) :: indexDotState
type(tDisloTungstenState), allocatable, dimension(:) :: state
2021-07-27 02:26:40 +05:30
type(tDisloTungstenDependentState), allocatable, dimension(:) :: dependentState
contains
2019-01-07 12:34:02 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Perform module initialization.
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
module function plastic_dislotungsten_init() result(myPlasticity)
2020-08-15 19:32:10 +05:30
logical, dimension(:), allocatable :: myPlasticity
2019-03-22 15:25:17 +05:30
integer :: &
ph, i, &
2021-03-05 01:46:36 +05:30
Nmembers, &
2019-03-22 15:25:17 +05:30
sizeState, sizeDotState, &
startIndex, endIndex
2020-03-16 19:32:40 +05:30
integer, dimension(:), allocatable :: &
2020-03-16 03:37:23 +05:30
N_sl
real(pREAL),dimension(:), allocatable :: &
f_edge, & !< edge character fraction of total dislocation density
2020-03-16 16:23:36 +05:30
rho_mob_0, & !< initial dislocation density
rho_dip_0 !< initial dipole density
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
myPlasticity = plastic_active('dislotungsten')
if (count(myPlasticity) == 0) return
2021-01-26 12:24:24 +05:30
print'(/,1x,a)', '<<<+- phase:mechanical:plastic:dislotungsten init -+>>>'
2021-02-13 17:37:35 +05:30
print'(/,1x,a)', 'D. Cereceda et al., International Journal of Plasticity 78:242256, 2016'
print'( 1x,a)', 'https://doi.org/10.1016/j.ijplas.2015.09.002'
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')
allocate(param(phases%length))
allocate(indexDotState(phases%length))
allocate(state(phases%length))
allocate(dependentState(phases%length))
extmsg = ''
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), dst => dependentState(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
2023-06-23 00:35:54 +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
2023-06-04 10:47:38 +05:30
prm%isotropic_bound = pl%get_asStr('isotropic_bound',defaultVal='isostrain')
2018-11-28 00:30:45 +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 19:32:40 +05:30
prm%sum_N_sl = sum(abs(N_sl))
2019-03-22 15:25:17 +05:30
slipActive: if (prm%sum_N_sl > 0) then
2023-07-15 23:58:57 +05:30
prm%systems_sl = crystal_labels_slip(N_sl,phase_lattice(ph))
prm%P_sl = 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
allocate(a_nS(3,size(pl%get_as1dReal('a_nonSchmid_110',defaultVal=emptyRealArray))))
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)
deallocate(a_nS)
2019-03-22 15:25:17 +05:30
else
prm%P_nS_pos = +prm%P_sl
prm%P_nS_neg = -prm%P_sl
2021-10-30 14:06:47 +05:30
end if
prm%dipoleformation = .not. pl%get_asBool('no_dipole_formation', defaultVal=.false.)
prm%D = pl%get_asReal('D')
prm%D_0 = pl%get_asReal('D_0')
prm%Q_cl = pl%get_asReal('Q_cl')
f_edge = math_expand(pl%get_as1dReal('f_edge', requiredSize=size(N_sl), &
defaultVal=[(0.5_pREAL, i=1,size(N_sl))]),N_sl)
rho_mob_0 = math_expand(pl%get_as1dReal('rho_mob_0', requiredSize=size(N_sl)),N_sl)
rho_dip_0 = math_expand(pl%get_as1dReal('rho_dip_0', requiredSize=size(N_sl)),N_sl)
prm%b_sl = math_expand(pl%get_as1dReal('b_sl', requiredSize=size(N_sl)),N_sl)
prm%Q_s = math_expand(pl%get_as1dReal('Q_s', requiredSize=size(N_sl)),N_sl)
prm%i_sl = math_expand(pl%get_as1dReal('i_sl', requiredSize=size(N_sl)),N_sl)
prm%tau_Peierls = math_expand(pl%get_as1dReal('tau_Peierls', requiredSize=size(N_sl)),N_sl)
prm%p = math_expand(pl%get_as1dReal('p_sl', requiredSize=size(N_sl)),N_sl)
prm%q = math_expand(pl%get_as1dReal('q_sl', requiredSize=size(N_sl)),N_sl)
prm%h = math_expand(pl%get_as1dReal('h', requiredSize=size(N_sl)),N_sl)
prm%w = math_expand(pl%get_as1dReal('w', requiredSize=size(N_sl)),N_sl)
prm%omega = math_expand(pl%get_as1dReal('omega', requiredSize=size(N_sl)),N_sl)
prm%B = math_expand(pl%get_as1dReal('B', requiredSize=size(N_sl)),N_sl)
prm%d_caron = prm%b_sl * pl%get_asReal('D_a')
prm%f_at = prm%b_sl**3*pl%get_asReal('f_at')
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))
prm%forestProjection = spread( f_edge,1,prm%sum_N_sl) &
2023-07-15 23:58:57 +05:30
* crystal_forestProjection_edge (N_sl,phase_lattice(ph),phase_cOverA(ph)) &
+ spread(1.0_pREAL-f_edge,1,prm%sum_N_sl) &
2023-07-15 23:58:57 +05:30
* crystal_forestProjection_screw(N_sl,phase_lattice(ph),phase_cOverA(ph))
2019-03-22 15:25:17 +05:30
! sanity checks
if ( prm%D_0 < 0.0_pREAL) extmsg = trim(extmsg)//' D_0'
if ( prm%Q_cl <= 0.0_pREAL) extmsg = trim(extmsg)//' Q_cl'
if (any(rho_mob_0 < 0.0_pREAL)) extmsg = trim(extmsg)//' rho_mob_0'
if (any(rho_dip_0 < 0.0_pREAL)) extmsg = trim(extmsg)//' rho_dip_0'
if (any(prm%b_sl <= 0.0_pREAL)) extmsg = trim(extmsg)//' b_sl'
if (any(prm%Q_s <= 0.0_pREAL)) extmsg = trim(extmsg)//' Q_s'
if (any(prm%tau_Peierls < 0.0_pREAL)) extmsg = trim(extmsg)//' tau_Peierls'
if (any(prm%B < 0.0_pREAL)) extmsg = trim(extmsg)//' B'
if (any(prm%d_caron < 0.0_pREAL)) extmsg = trim(extmsg)//' d_caron(D_a,b_sl)'
if (any(prm%f_at <= 0.0_pREAL)) extmsg = trim(extmsg)//' f_at or b_sl'
2019-03-22 15:25:17 +05:30
else slipActive
rho_mob_0 = emptyRealArray
rho_dip_0 = emptyRealArray
allocate(prm%b_sl, &
prm%d_caron, &
prm%i_sl, &
prm%f_at, &
prm%tau_Peierls, &
prm%Q_s, &
prm%p, &
prm%q, &
prm%B, &
prm%h, &
prm%w, &
prm%omega, &
2020-03-16 16:23:36 +05:30
source = emptyRealArray)
allocate(prm%forestProjection(0,0))
allocate(prm%h_sl_sl (0,0))
2021-10-30 14:06:47 +05:30
end if slipActive
2018-11-28 10:29:03 +05:30
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2023-01-23 13:01:59 +05:30
Nmembers = count(material_ID_phase == ph)
2019-03-22 15:25:17 +05:30
sizeDotState = size(['rho_mob ','rho_dip ','gamma_sl']) * prm%sum_N_sl
sizeState = sizeDotState
2021-03-05 01:46:36 +05:30
call phase_allocateState(plasticState(ph),Nmembers,sizeState,sizeDotState,0)
deallocate(plasticState(ph)%dotState) ! ToDo: remove dotState completely
2018-12-05 03:00:07 +05:30
!--------------------------------------------------------------------------------------------------
2020-03-16 03:37:23 +05:30
! state aliases and initialization
2019-03-22 15:25:17 +05:30
startIndex = 1
endIndex = prm%sum_N_sl
idx_dot%rho_mob = [startIndex,endIndex]
2021-02-16 20:36:09 +05:30
stt%rho_mob => plasticState(ph)%state(startIndex:endIndex,:)
stt%rho_mob = spread(rho_mob_0,2,Nmembers)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_rho',defaultVal=1.0_pREAL)
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_rho'
2019-03-22 15:25:17 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_sl
idx_dot%rho_dip = [startIndex,endIndex]
2021-02-16 20:36:09 +05:30
stt%rho_dip => plasticState(ph)%state(startIndex:endIndex,:)
stt%rho_dip = spread(rho_dip_0,2,Nmembers)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_rho',defaultVal=1.0_pREAL)
2019-03-22 15:25:17 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_sl
idx_dot%gamma_sl = [startIndex,endIndex]
2021-02-16 20:36:09 +05:30
stt%gamma_sl => 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'
allocate(dst%Lambda_sl(prm%sum_N_sl,Nmembers), source=0.0_pREAL)
allocate(dst%tau_pass (prm%sum_N_sl,Nmembers), source=0.0_pREAL)
2019-03-22 15:25:17 +05:30
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
2021-10-30 14:06:47 +05:30
end do
2018-12-04 04:36:46 +05:30
end function plastic_dislotungsten_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
2021-01-26 12:03:04 +05:30
pure module subroutine dislotungsten_LpAndItsTangent(Lp,dLp_dMp, &
Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2019-03-22 15:25:17 +05:30
Lp !< plastic velocity gradient
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2019-03-22 15:25:17 +05:30
dLp_dMp !< derivative of Lp with respect to the Mandel stress
real(pREAL), dimension(3,3), intent(in) :: &
2019-03-22 15:25:17 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2019-03-22 15:25:17 +05:30
integer :: &
i,k,l,m,n
real(pREAL) :: &
T !< temperature
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
dot_gamma, ddot_gamma_dtau
real(pREAL), dimension(3,3,param(ph)%sum_N_sl) :: &
P_nS
2021-10-30 14:06:47 +05:30
T = thermal_T(ph,en)
Lp = 0.0_pREAL
dLp_dMp = 0.0_pREAL
associate(prm => param(ph))
call kinetics(Mp,T,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
2021-10-30 14:06:47 +05:30
do i = 1, prm%sum_N_sl
Lp = Lp + dot_gamma(i)*prm%P_sl(1:3,1:3,i)
2021-10-30 14:06:47 +05:30
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_sl(k,l,i) * P_nS(m,n,i)
2021-10-30 14:06:47 +05:30
end do
2019-03-22 15:25:17 +05:30
end associate
2021-01-26 12:03:04 +05:30
end subroutine dislotungsten_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief Calculate the rate of change of microstructure.
!--------------------------------------------------------------------------------------------------
module function dislotungsten_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2020-03-16 04:42:18 +05:30
Mp !< Mandel stress
2019-03-22 15:25:17 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
tau_eff, &
2019-03-22 15:36:08 +05:30
v_cl, &
dot_rho_dip_formation, &
dot_rho_dip_climb, &
2021-07-24 01:34:39 +05:30
d_hat
real(pREAL) :: &
mu, nu, T
2021-10-30 14:06:47 +05:30
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph), &
dot_rho_mob => dotState(indexDotState(ph)%rho_mob(1):indexDotState(ph)%rho_mob(2)), &
dot_rho_dip => dotState(indexDotState(ph)%rho_dip(1):indexDotState(ph)%rho_dip(2)), &
dot_gamma => dotState(indexDotState(ph)%gamma_sl(1):indexDotState(ph)%gamma_sl(2)))
mu = elastic_mu(ph,en,prm%isotropic_bound)
nu = elastic_nu(ph,en,prm%isotropic_bound)
T = thermal_T(ph,en)
2021-10-30 14:06:47 +05:30
call kinetics(Mp,T,ph,en,&
dot_gamma, tau = tau_eff)
2021-10-30 14:06:47 +05:30
dot_gamma = abs(dot_gamma)
2021-10-30 14:06:47 +05:30
where(dEq0(dot_gamma))
dot_rho_dip_formation = 0.0_pREAL
dot_rho_dip_climb = 0.0_pREAL
2021-10-30 14:06:47 +05:30
else where
d_hat = math_clip(mu*prm%b_sl/(8.0_pREAL*PI*(1.0_pREAL-nu)*tau_eff), &
2023-09-16 07:11:37 +05:30
left = prm%d_caron, & ! lower limit
right = dst%Lambda_sl(:,en)) ! upper limit
dot_rho_dip_formation = merge(dot_gamma * 2.0_pREAL*(d_hat-prm%d_caron)/prm%b_sl * stt%rho_mob(:,en), &
0.0_pREAL, &
2021-10-30 14:06:47 +05:30
prm%dipoleformation)
v_cl = (3.0_pREAL*mu*prm%D_0*exp(-prm%Q_cl/(K_B*T))*prm%f_at/(2.0_pREAL*PI*K_B*T)) &
* (1.0_pREAL/(d_hat+prm%d_caron))
dot_rho_dip_climb = (4.0_pREAL*v_cl*stt%rho_dip(:,en))/(d_hat-prm%d_caron) ! ToDo: Discuss with Franz: Stress dependency?
2021-10-30 14:06:47 +05:30
end where
dot_rho_mob = dot_gamma / (prm%b_sl*dst%Lambda_sl(:,en)) & ! multiplication
- dot_rho_dip_formation &
- dot_gamma * 2.0_pREAL*prm%d_caron/prm%b_sl * stt%rho_mob(:,en) ! spontaneous annihilation of 2 edges
dot_rho_dip = dot_rho_dip_formation &
- dot_rho_dip_climb &
- dot_gamma * 2.0_pREAL*prm%d_caron/prm%b_sl * stt%rho_dip(:,en) ! spontaneous annihilation of an edge with a dipole
2019-03-22 15:25:17 +05:30
end associate
2018-12-05 03:00:07 +05:30
end function dislotungsten_dotState
2018-12-05 03:00:07 +05:30
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate derived quantities from state.
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine dislotungsten_dependentState(ph,en)
2021-10-30 14:06:47 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
Lambda_sl_inv
2021-11-18 17:16:37 +05:30
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
dst%tau_pass(:,en) = elastic_mu(ph,en,prm%isotropic_bound)*prm%b_sl &
2021-07-24 01:26:34 +05:30
* sqrt(matmul(prm%h_sl_sl,stt%rho_mob(:,en)+stt%rho_dip(:,en)))
Lambda_sl_inv = 1.0_pREAL/prm%D &
+ sqrt(matmul(prm%forestProjection,stt%rho_mob(:,en)+stt%rho_dip(:,en)))/prm%i_sl
dst%Lambda_sl(:,en) = Lambda_sl_inv**(-1.0_pREAL)
2019-03-22 15:25:17 +05:30
end associate
2019-01-07 12:06:11 +05:30
2021-01-26 16:47:00 +05:30
end subroutine dislotungsten_dependentState
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
module subroutine plastic_dislotungsten_result(ph,group)
integer, intent(in) :: ph
2019-04-04 11:22:36 +05:30
character(len=*), intent(in) :: group
integer :: ou
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
do ou = 1,size(prm%output)
select case(trim(prm%output(ou)))
case('rho_mob')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%rho_mob,group,trim(prm%output(ou)), &
'mobile dislocation density','1/m²',prm%systems_sl)
case('rho_dip')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%rho_dip,group,trim(prm%output(ou)), &
'dislocation dipole density','1/m²',prm%systems_sl)
case('gamma_sl')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%gamma_sl,group,trim(prm%output(ou)), &
'plastic shear','1',prm%systems_sl)
case('Lambda_sl')
2023-01-19 22:07:45 +05:30
call result_writeDataset(dst%Lambda_sl,group,trim(prm%output(ou)), &
'mean free path for slip','m',prm%systems_sl)
case('tau_pass')
2023-01-19 22:07:45 +05:30
call result_writeDataset(dst%tau_pass,group,trim(prm%output(ou)), &
'threshold stress for slip','Pa',prm%systems_sl)
end select
2021-10-30 14:06:47 +05:30
end do
end associate
2023-01-19 22:07:45 +05:30
end subroutine plastic_dislotungsten_result
2019-01-07 11:54:02 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate shear rates on slip systems, their derivatives with respect to resolved
! stress, and the resolved stress.
2019-01-07 11:54:02 +05:30
!> @details Derivatives and resolved stress 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.
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
pure subroutine kinetics(Mp,T,ph,en, &
dot_gamma,ddot_gamma_dtau,tau)
real(pREAL), dimension(3,3), intent(in) :: &
2019-03-22 15:25:17 +05:30
Mp !< Mandel stress
real(pREAL), intent(in) :: &
2019-03-22 15:25:17 +05:30
T !< temperature
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
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, &
tau
2021-10-30 14:06:47 +05:30
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
StressRatio, &
StressRatio_p,StressRatio_pminus1, &
2021-10-31 02:32:51 +05:30
tau_pos, tau_neg, tau_eff, &
t_n,t_k, dtk,dtn
integer :: i
2021-10-30 14:06:47 +05:30
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
tau_pos = [(math_tensordot(Mp,prm%P_nS_pos(1:3,1:3,i)),i=1,prm%sum_N_sl)]
tau_neg = [(math_tensordot(Mp,prm%P_nS_neg(1:3,1:3,i)),i=1,prm%sum_N_sl)]
2023-09-22 02:02:38 +05:30
tau_eff = math_clip(max(tau_pos,tau_neg) - dst%tau_pass(:,en),left = 0.0_pREAL)
if (present(tau)) tau = tau_eff
associate(BoltzmannRatio => prm%Q_s/(K_B*T), &
b_rho => stt%rho_mob(:,en) * prm%b_sl, &
2021-10-30 14:06:47 +05:30
effectiveLength => dst%Lambda_sl(:,en) - prm%w)
2021-10-31 02:32:51 +05:30
where(tau_eff > tol_math_check)
2021-10-31 02:32:51 +05:30
StressRatio = tau_eff/prm%tau_Peierls
2021-10-30 14:06:47 +05:30
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pREAL)
t_n = prm%b_sl*exp(BoltzmannRatio*(1.0_pREAL-StressRatio_p) ** prm%q) &
2021-10-31 01:38:55 +05:30
/ (prm%omega*effectiveLength)
t_k = effectiveLength * prm%B /(2.0_pREAL*prm%b_sl*tau_eff) ! corrected eq. (14)
dot_gamma = b_rho * prm%h/(t_n + t_k) * merge(+1.0_pREAL,-1.0_pREAL, tau_pos>tau_neg)
else where
dot_gamma = 0.0_pREAL
end where
if (present(ddot_gamma_dtau)) then
where(tau_eff > tol_math_check)
dtn = -1.0_pREAL * t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pREAL-StressRatio_p)**(prm%q - 1.0_pREAL) &
2021-10-31 21:59:35 +05:30
* StressRatio_pminus1 / prm%tau_Peierls
dtk = -1.0_pREAL * t_k / tau_eff
ddot_gamma_dtau = -1.0_pREAL * dot_gamma * (dtn + dtk) / (t_n + t_k)
else where
ddot_gamma_dtau = 0.0_pREAL
end where
2021-10-30 14:06:47 +05:30
end if
2021-10-30 14:06:47 +05:30
end associate
2019-03-22 15:25:17 +05:30
end associate
2018-12-05 02:03:32 +05:30
end subroutine kinetics
2021-01-26 05:41:32 +05:30
end submodule dislotungsten