DAMASK_EICMD/src/phase_mechanical_plastic_di...

537 lines
25 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
real(pReal), parameter :: &
2019-03-22 15:25:17 +05:30
kB = 1.38e-23_pReal !< Boltzmann constant in J/Kelvin
2019-12-03 02:21:25 +05:30
type :: tParameters
2019-03-22 15:25:17 +05:30
real(pReal) :: &
2020-03-16 16:23:36 +05:30
D = 1.0_pReal, & !< grain size
mu = 1.0_pReal, & !< equivalent shear modulus
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]
v_0, & !< dislocation velocity prefactor [m/s]
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
character(len=pStringLen), allocatable, dimension(:) :: &
output
2019-03-22 15:25:17 +05:30
logical :: &
dipoleFormation !< flag indicating consideration of dipole formation
end type !< container type for internal constitutive parameters
2020-08-15 19:32:10 +05:30
type :: tDisloTungstenState
2019-03-22 15:25:17 +05:30
real(pReal), dimension(:,:), pointer :: &
rho_mob, &
rho_dip, &
gamma_sl
2020-08-15 19:32:10 +05:30
end type tDisloTungstenState
2020-08-15 19:32:10 +05:30
type :: tDisloTungstendependentState
2019-03-22 15:25:17 +05:30
real(pReal), dimension(:,:), allocatable :: &
Lambda_sl, &
2021-07-24 01:26:34 +05:30
tau_pass
2020-08-15 19:32:10 +05:30
end type tDisloTungstendependentState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
2019-12-03 02:21:25 +05:30
type(tParameters), allocatable, dimension(:) :: param
2020-08-15 19:32:10 +05:30
type(tDisloTungstenState), allocatable, dimension(:) :: &
2019-03-22 15:25:17 +05:30
dotState, &
state
2020-08-15 19:32:10 +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 :: &
2021-02-16 20:36:09 +05:30
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
2020-03-16 19:32:40 +05:30
real(pReal),dimension(:), allocatable :: &
2020-03-16 16:23:36 +05:30
rho_mob_0, & !< initial dislocation density
rho_dip_0, & !< initial dipole density
a !< non-Schmid coefficients
2019-03-22 15:25:17 +05:30
character(len=pStringLen) :: &
extmsg = ''
2020-08-15 19:32:10 +05:30
class(tNode), pointer :: &
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
2021-02-13 23:27:41 +05:30
print'(/,a)', ' <<<+- phase:mechanical:plastic:dislotungsten init -+>>>'
print'(a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
2021-02-13 17:37:35 +05:30
2021-03-18 20:06:40 +05:30
print*, 'D. Cereceda et al., International Journal of Plasticity 78:242256, 2016'
2021-03-26 01:33:46 +05:30
print*, 'https://doi.org/10.1016/j.ijplas.2015.09.002'
phases => config_material%get('phase')
allocate(param(phases%length))
allocate(state(phases%length))
allocate(dotState(phases%length))
allocate(dependentState(phases%length))
2021-02-16 20:36:09 +05:30
do ph = 1, phases%length
if(.not. myPlasticity(ph)) cycle
associate(prm => param(ph), dot => dotState(ph), stt => state(ph), dst => dependentState(ph))
phase => phases%get(ph)
2021-03-25 23:52:59 +05:30
mech => phase%get('mechanical')
pl => mech%get('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
prm%mu = elastic_mu(ph)
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
2021-07-21 19:16:38 +05:30
prm%P_sl = 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)
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)
2019-03-22 15:25:17 +05:30
else
2021-07-21 19:59:57 +05:30
prm%P_nS_pos = prm%P_sl
prm%P_nS_neg = prm%P_sl
2019-03-22 15:25:17 +05:30
endif
2021-06-19 18:43:27 +05:30
prm%h_sl_sl = lattice_interaction_SlipBySlip(N_sl,pl%get_as1dFloat('h_sl-sl'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
prm%forestProjection = lattice_forestProjection_edge(N_sl,phase_lattice(ph),&
phase_cOverA(ph))
2020-03-16 16:23:36 +05:30
prm%forestProjection = transpose(prm%forestProjection)
2021-03-11 22:30:07 +05:30
rho_mob_0 = pl%get_as1dFloat('rho_mob_0', requiredSize=size(N_sl))
rho_dip_0 = pl%get_as1dFloat('rho_dip_0', requiredSize=size(N_sl))
prm%v_0 = pl%get_as1dFloat('v_0', requiredSize=size(N_sl))
prm%b_sl = pl%get_as1dFloat('b_sl', requiredSize=size(N_sl))
prm%Q_s = pl%get_as1dFloat('Q_s', requiredSize=size(N_sl))
2020-03-16 03:37:23 +05:30
2021-03-11 22:30:07 +05:30
prm%i_sl = pl%get_as1dFloat('i_sl', requiredSize=size(N_sl))
prm%tau_Peierls = pl%get_as1dFloat('tau_Peierls', requiredSize=size(N_sl))
2021-07-23 10:52:52 +05:30
prm%p = pl%get_as1dFloat('p_sl', requiredSize=size(N_sl))
prm%q = pl%get_as1dFloat('q_sl', requiredSize=size(N_sl))
2021-03-11 22:30:07 +05:30
prm%h = pl%get_as1dFloat('h', requiredSize=size(N_sl))
prm%w = pl%get_as1dFloat('w', requiredSize=size(N_sl))
prm%omega = pl%get_as1dFloat('omega', requiredSize=size(N_sl))
prm%B = pl%get_as1dFloat('B', requiredSize=size(N_sl))
2020-08-15 19:32:10 +05:30
2021-07-23 10:52:52 +05:30
prm%D = pl%get_asFloat('D')
prm%D_0 = pl%get_asFloat('D_0')
prm%Q_cl = pl%get_asFloat('Q_cl')
prm%f_at = pl%get_asFloat('f_at') * prm%b_sl**3.0_pReal
prm%dipoleformation = .not. pl%get_asBool('no_dipole_formation', defaultVal = .false.)
2019-03-22 15:25:17 +05:30
! expand: family => system
2021-07-23 10:52:52 +05:30
rho_mob_0 = math_expand(rho_mob_0, N_sl)
rho_dip_0 = math_expand(rho_dip_0, N_sl)
prm%q = math_expand(prm%q, N_sl)
prm%p = math_expand(prm%p, N_sl)
prm%Q_s = math_expand(prm%Q_s, N_sl)
prm%b_sl = math_expand(prm%b_sl, N_sl)
prm%h = math_expand(prm%h, N_sl)
prm%w = math_expand(prm%w, N_sl)
prm%omega = math_expand(prm%omega, N_sl)
prm%tau_Peierls = math_expand(prm%tau_Peierls, N_sl)
prm%v_0 = math_expand(prm%v_0, N_sl)
prm%B = math_expand(prm%B, N_sl)
prm%i_sl = math_expand(prm%i_sl, N_sl)
prm%f_at = math_expand(prm%f_at, N_sl)
2021-07-24 01:34:39 +05:30
prm%d_caron = pl%get_asFloat('D_a') * prm%b_sl
2019-03-22 15:25:17 +05:30
! sanity checks
2020-03-16 16:23:36 +05:30
if ( prm%D_0 <= 0.0_pReal) extmsg = trim(extmsg)//' D_0'
if ( prm%Q_cl <= 0.0_pReal) extmsg = trim(extmsg)//' Q_cl'
2020-08-15 19:32:10 +05:30
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'
2020-09-22 18:05:05 +05:30
if (any(prm%v_0 < 0.0_pReal)) extmsg = trim(extmsg)//' v_0'
2020-03-16 16:23:36 +05:30
if (any(prm%b_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' b_sl'
2020-09-22 18:05:05 +05:30
if (any(prm%Q_s <= 0.0_pReal)) extmsg = trim(extmsg)//' Q_s'
2020-09-23 04:25:19 +05:30
if (any(prm%tau_Peierls < 0.0_pReal)) extmsg = trim(extmsg)//' tau_Peierls'
2021-07-24 01:34:39 +05:30
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)'
2020-09-22 18:05:05 +05:30
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
2020-03-16 16:23:36 +05:30
rho_mob_0= emptyRealArray; rho_dip_0 = emptyRealArray
2021-07-24 01:34:39 +05:30
allocate(prm%b_sl,prm%d_caron,prm%i_sl,prm%f_at,prm%tau_Peierls, &
2020-09-22 18:05:05 +05:30
prm%Q_s,prm%v_0,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))
2019-03-22 15:25:17 +05:30
endif slipActive
2018-11-28 10:29:03 +05:30
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2021-04-06 15:08:44 +05:30
Nmembers = count(material_phaseID == 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)
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
2021-02-16 20:36:09 +05:30
stt%rho_mob => plasticState(ph)%state(startIndex:endIndex,:)
2021-03-05 01:46:36 +05:30
stt%rho_mob = spread(rho_mob_0,2,Nmembers)
2021-02-16 20:36:09 +05:30
dot%rho_mob => plasticState(ph)%dotState(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('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
2021-02-16 20:36:09 +05:30
stt%rho_dip => plasticState(ph)%state(startIndex:endIndex,:)
2021-03-05 01:46:36 +05:30
stt%rho_dip = spread(rho_dip_0,2,Nmembers)
2021-02-16 20:36:09 +05:30
dot%rho_dip => plasticState(ph)%dotState(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_rho',defaultVal=1.0_pReal)
2019-03-22 15:25:17 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_sl
2021-02-16 20:36:09 +05:30
stt%gamma_sl => plasticState(ph)%state(startIndex:endIndex,:)
dot%gamma_sl => plasticState(ph)%dotState(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_gamma',defaultVal=1.0e-6_pReal)
if(any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_gamma'
2021-07-24 01:26:34 +05:30
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)//'(dislotungsten)')
2020-03-15 00:33:57 +05:30
2019-03-22 15:25:17 +05:30
enddo
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, &
2021-04-25 11:36:52 +05:30
Mp,T,ph,en)
2019-03-22 15:25:17 +05:30
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
2019-03-22 15:25:17 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
T !< temperature
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), dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
dot_gamma_pos,dot_gamma_neg, &
ddot_gamma_dtau_pos,ddot_gamma_dtau_neg
2019-03-22 15:25:17 +05:30
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
associate(prm => param(ph))
2021-04-25 11:36:52 +05:30
call kinetics(Mp,T,ph,en,dot_gamma_pos,dot_gamma_neg,ddot_gamma_dtau_pos,ddot_gamma_dtau_neg)
2019-03-22 15:25:17 +05:30
do i = 1, prm%sum_N_sl
2020-03-16 18:03:14 +05:30
Lp = Lp + (dot_gamma_pos(i)+dot_gamma_neg(i))*prm%P_sl(1:3,1:3,i)
2019-03-22 15:25:17 +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) &
2021-07-21 19:59:57 +05:30
+ ddot_gamma_dtau_pos(i) * prm%P_sl(k,l,i) * prm%P_nS_pos(m,n,i) &
+ ddot_gamma_dtau_neg(i) * prm%P_sl(k,l,i) * prm%P_nS_neg(m,n,i)
2019-03-22 15:25:17 +05:30
enddo
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.
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine dislotungsten_dotState(Mp,T,ph,en)
2019-03-22 15:25:17 +05:30
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
real(pReal), intent(in) :: &
2020-03-16 04:42:18 +05:30
T !< temperature
2019-03-22 15:25:17 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_pos, dot_gamma_neg,&
2019-03-22 15:25:17 +05:30
tau_pos,&
tau_neg, &
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
2021-07-24 01:34:39 +05:30
associate(prm => param(ph), stt => state(ph), dot => dotState(ph), dst => dependentState(ph))
2021-04-25 11:36:52 +05:30
call kinetics(Mp,T,ph,en,&
2021-07-21 20:04:33 +05:30
dot_gamma_pos,dot_gamma_neg, &
2019-03-22 15:25:17 +05:30
tau_pos_out = tau_pos,tau_neg_out = tau_neg)
2021-07-21 20:04:33 +05:30
dot%gamma_sl(:,en) = (dot_gamma_pos+dot_gamma_neg) ! ToDo: needs to be abs
2019-03-22 15:25:17 +05:30
where(dEq0(tau_pos)) ! ToDo: use avg of pos and neg
2019-03-22 15:36:08 +05:30
dot_rho_dip_formation = 0.0_pReal
dot_rho_dip_climb = 0.0_pReal
2019-03-22 15:25:17 +05:30
else where
2021-07-24 01:34:39 +05:30
d_hat = math_clip(3.0_pReal*prm%mu*prm%b_sl/(16.0_pReal*PI*abs(tau_pos)), &
prm%d_caron, & ! lower limit
dst%Lambda_sl(:,en)) ! upper limit
dot_rho_dip_formation = merge(2.0_pReal*d_hat* stt%rho_mob(:,en)*abs(dot%gamma_sl(:,en))/prm%b_sl, & ! ToDo: ignore region of spontaneous annihilation
2019-03-22 15:36:08 +05:30
0.0_pReal, &
prm%dipoleformation)
2021-07-24 01:34:39 +05:30
v_cl = (3.0_pReal*prm%mu*prm%D_0*exp(-prm%Q_cl/(kB*T))*prm%f_at/(2.0_pReal*PI*kB*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?
2019-03-22 15:25:17 +05:30
end where
2021-04-25 11:36:52 +05:30
dot%rho_mob(:,en) = abs(dot%gamma_sl(:,en))/(prm%b_sl*dst%Lambda_sl(:,en)) & ! multiplication
2019-03-22 15:36:08 +05:30
- dot_rho_dip_formation &
2021-07-24 01:34:39 +05:30
- (2.0_pReal*prm%d_caron)/prm%b_sl*stt%rho_mob(:,en)*abs(dot%gamma_sl(:,en)) ! Spontaneous annihilation of 2 edges
2021-04-25 11:36:52 +05:30
dot%rho_dip(:,en) = dot_rho_dip_formation &
2021-07-24 01:34:39 +05:30
- (2.0_pReal*prm%d_caron)/prm%b_sl*stt%rho_dip(:,en)*abs(dot%gamma_sl(:,en)) & ! Spontaneous annihilation of an edge with a dipole
2019-03-22 15:36:08 +05:30
- dot_rho_dip_climb
2019-03-22 15:25:17 +05:30
end associate
2018-12-05 03:00:07 +05:30
2021-01-26 16:47:00 +05:30
end subroutine 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)
2019-03-22 15:25:17 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
dislocationSpacing
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph),dst => dependentState(ph))
2021-07-23 03:39:51 +05:30
dislocationSpacing = sqrt(matmul(prm%forestProjection,stt%rho_mob(:,en)+stt%rho_dip(:,en)))
2021-07-24 01:26:34 +05:30
dst%tau_pass(:,en) = prm%mu*prm%b_sl &
* sqrt(matmul(prm%h_sl_sl,stt%rho_mob(:,en)+stt%rho_dip(:,en)))
2021-07-23 03:39:51 +05:30
dst%Lambda_sl(:,en) = prm%D/(1.0_pReal+prm%D*dislocationSpacing/prm%i_sl)
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.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_dislotungsten_results(ph,group)
integer, intent(in) :: ph
2019-04-04 11:22:36 +05:30
character(len=*), intent(in) :: group
integer :: o
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
case('rho_mob')
2021-06-01 20:35:13 +05:30
if(prm%sum_N_sl>0) call results_writeDataset(stt%rho_mob,group,trim(prm%output(o)), &
'mobile dislocation density','1/m²')
case('rho_dip')
2021-06-01 20:35:13 +05:30
if(prm%sum_N_sl>0) call results_writeDataset(stt%rho_dip,group,trim(prm%output(o)), &
2021-05-20 21:28:42 +05:30
'dislocation dipole density','1/m²')
case('gamma_sl')
2021-06-01 20:35:13 +05:30
if(prm%sum_N_sl>0) call results_writeDataset(stt%gamma_sl,group,trim(prm%output(o)), &
'plastic shear','1')
case('Lambda_sl')
2021-06-01 20:35:13 +05:30
if(prm%sum_N_sl>0) call results_writeDataset(dst%Lambda_sl,group,trim(prm%output(o)), &
'mean free path for slip','m')
case('tau_pass')
2021-07-24 01:26:34 +05:30
if(prm%sum_N_sl>0) call results_writeDataset(dst%tau_pass,group,trim(prm%output(o)), &
'threshold stress for slip','Pa')
end select
enddo outputsLoop
end associate
end subroutine plastic_dislotungsten_results
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.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
pure subroutine kinetics(Mp,T,ph,en, &
2019-03-22 15:19:24 +05:30
dot_gamma_pos,dot_gamma_neg,ddot_gamma_dtau_pos,ddot_gamma_dtau_neg,tau_pos_out,tau_neg_out)
2019-03-22 15:25:17 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
T !< temperature
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pReal), intent(out), dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
dot_gamma_pos, &
dot_gamma_neg
real(pReal), intent(out), optional, dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
ddot_gamma_dtau_pos, &
ddot_gamma_dtau_neg, &
tau_pos_out, &
tau_neg_out
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2019-03-22 15:25:17 +05:30
StressRatio, &
StressRatio_p,StressRatio_pminus1, &
dvel, vel, &
tau_pos,tau_neg, &
2021-07-24 01:26:34 +05:30
t_n, t_k, dtk,dtn
2019-03-22 15:25:17 +05:30
integer :: j
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
2019-03-22 15:25:17 +05:30
do j = 1, prm%sum_N_sl
2021-07-21 19:59:57 +05:30
tau_pos(j) = math_tensordot(Mp,prm%P_nS_pos(1:3,1:3,j))
tau_neg(j) = math_tensordot(Mp,prm%P_nS_neg(1:3,1:3,j))
2019-03-22 15:25:17 +05:30
enddo
2019-03-22 15:25:17 +05:30
if (present(tau_pos_out)) tau_pos_out = tau_pos
if (present(tau_neg_out)) tau_neg_out = tau_neg
2020-09-22 18:05:05 +05:30
associate(BoltzmannRatio => prm%Q_s/(kB*T), &
2021-04-25 11:36:52 +05:30
dot_gamma_0 => stt%rho_mob(:,en)*prm%b_sl*prm%v_0, &
effectiveLength => dst%Lambda_sl(:,en) - prm%w)
2021-07-24 01:26:34 +05:30
significantPositiveTau: where(abs(tau_pos)-dst%tau_pass(:,en) > tol_math_check)
StressRatio = (abs(tau_pos)-dst%tau_pass(:,en))/prm%tau_Peierls
2019-03-22 15:25:17 +05:30
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
2021-07-24 01:26:34 +05:30
t_n = prm%b_sl/(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)*prm%omega*effectiveLength)
2019-09-26 00:51:57 +05:30
t_k = effectiveLength * prm%B /(2.0_pReal*prm%b_sl*tau_pos)
2020-09-22 18:05:05 +05:30
vel = prm%h/(t_n + t_k)
2019-03-22 15:25:17 +05:30
dot_gamma_pos = dot_gamma_0 * sign(vel,tau_pos) * 0.5_pReal
else where significantPositiveTau
dot_gamma_pos = 0.0_pReal
end where significantPositiveTau
2019-03-22 15:25:17 +05:30
if (present(ddot_gamma_dtau_pos)) then
2021-07-24 01:26:34 +05:30
significantPositiveTau2: where(abs(tau_pos)-dst%tau_pass(:,en) > tol_math_check)
2019-09-26 00:51:57 +05:30
dtn = -1.0_pReal * t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
2020-09-23 04:25:19 +05:30
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_Peierls
2019-09-26 00:51:57 +05:30
dtk = -1.0_pReal * t_k / tau_pos
2020-09-22 18:05:05 +05:30
dvel = -1.0_pReal * prm%h * (dtk + dtn) / (t_n + t_k)**2.0_pReal
2019-09-26 00:51:57 +05:30
ddot_gamma_dtau_pos = dot_gamma_0 * dvel* 0.5_pReal
else where significantPositiveTau2
ddot_gamma_dtau_pos = 0.0_pReal
end where significantPositiveTau2
2019-03-22 15:25:17 +05:30
endif
2021-07-24 01:26:34 +05:30
significantNegativeTau: where(abs(tau_neg)-dst%tau_pass(:,en) > tol_math_check)
StressRatio = (abs(tau_neg)-dst%tau_pass(:,en))/prm%tau_Peierls
2019-03-22 15:25:17 +05:30
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
2021-07-24 01:26:34 +05:30
t_n = prm%b_sl/(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)*prm%omega*effectiveLength)
2019-09-26 00:51:57 +05:30
t_k = effectiveLength * prm%B /(2.0_pReal*prm%b_sl*tau_pos)
2020-09-22 18:05:05 +05:30
vel = prm%h/(t_n + t_k)
2019-03-22 15:25:17 +05:30
dot_gamma_neg = dot_gamma_0 * sign(vel,tau_neg) * 0.5_pReal
else where significantNegativeTau
dot_gamma_neg = 0.0_pReal
end where significantNegativeTau
2019-03-22 15:25:17 +05:30
if (present(ddot_gamma_dtau_neg)) then
2021-07-24 01:26:34 +05:30
significantNegativeTau2: where(abs(tau_neg)-dst%tau_pass(:,en) > tol_math_check)
2019-09-26 00:51:57 +05:30
dtn = -1.0_pReal * t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
2020-09-23 04:25:19 +05:30
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_Peierls
2019-09-26 00:51:57 +05:30
dtk = -1.0_pReal * t_k / tau_neg
2020-09-22 18:05:05 +05:30
dvel = -1.0_pReal * prm%h * (dtk + dtn) / (t_n + t_k)**2.0_pReal
2019-09-26 00:51:57 +05:30
ddot_gamma_dtau_neg = dot_gamma_0 * dvel * 0.5_pReal
else where significantNegativeTau2
ddot_gamma_dtau_neg = 0.0_pReal
end where significantNegativeTau2
2019-03-22 15:25:17 +05:30
end if
2019-03-22 15:25:17 +05:30
end associate
end associate
2018-12-05 02:03:32 +05:30
end subroutine kinetics
2021-01-26 05:41:32 +05:30
end submodule dislotungsten