2015-04-21 20:46:13 +05:30
!--------------------------------------------------------------------------------------------------
2014-11-06 00:41:09 +05:30
!> @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
2019-01-06 04:03:18 +05:30
!> @brief crystal plasticity model for bcc metals, especially Tungsten
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2021-01-27 01:22:48 +05:30
submodule ( phase : plastic ) dislotungsten
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
real ( pReal ) , allocatable , dimension ( : ) :: &
2021-01-26 13:09:17 +05:30
b_sl , & !< magnitude me Burgers vector [m]
2019-03-22 15:25:17 +05:30
D_a , &
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
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
real ( pReal ) , allocatable , dimension ( : , : , : ) :: &
2020-03-16 18:03:14 +05:30
P_sl , &
2019-03-22 15:25:17 +05:30
nonSchmid_pos , &
nonSchmid_neg
integer :: &
sum_N_sl !< total number of active slip system
2020-02-14 13:30:14 +05:30
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-02-14 13:30:14 +05:30
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-02-14 13:30:14 +05:30
2020-08-15 19:32:10 +05:30
type :: tDisloTungstendependentState
2019-03-22 15:25:17 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable :: &
Lambda_sl , &
threshold_stress
2020-08-15 19:32:10 +05:30
end type tDisloTungstendependentState
2018-11-29 12:44:20 +05:30
2019-01-13 21:33:49 +05:30
!--------------------------------------------------------------------------------------------------
! 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
2018-11-28 21:42:06 +05:30
2014-11-06 00:41:09 +05:30
contains
2019-01-07 12:34:02 +05:30
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +05:30
!> @brief Perform module initialization.
2014-11-06 00:41:09 +05:30
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2020-11-29 16:28:39 +05:30
module function plastic_dislotungsten_init ( ) result ( myPlasticity )
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
2020-11-29 16:28:39 +05:30
myPlasticity = plastic_active ( 'dislotungsten' )
2021-02-14 11:36:14 +05:30
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 -+>>>'
2021-02-14 11:36:14 +05:30
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:242– 256, 2016'
2020-09-14 00:30:34 +05:30
print * , 'https://dx.doi.org/10.1016/j.ijplas.2015.09.002'
2020-02-14 13:30:14 +05:30
2020-09-13 14:09:17 +05:30
phases = > config_material % get ( 'phase' )
2021-02-14 11:36:14 +05:30
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__)
prm % output = output_asStrings ( pl )
#else
prm % output = pl % get_asStrings ( 'output' , defaultVal = emptyStringArray )
#endif
2020-02-14 13:30:14 +05:30
2020-03-15 00:33:57 +05:30
! This data is read in already in lattice
2021-02-16 20:36:09 +05:30
prm % mu = lattice_mu ( ph )
2019-01-06 04:03:18 +05:30
2018-11-28 00:30:45 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
2020-08-15 19:32:10 +05:30
N_sl = pl % get_asInts ( '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
2020-08-15 19:32:10 +05:30
prm % P_sl = lattice_SchmidMatrix_slip ( N_sl , phase % get_asString ( 'lattice' ) , &
phase % get_asFloat ( 'c/a' , defaultVal = 0.0_pReal ) )
2020-02-14 13:30:14 +05:30
2020-11-29 16:28:39 +05:30
if ( trim ( phase % get_asString ( 'lattice' ) ) == 'cI' ) then
2020-09-23 04:25:19 +05:30
a = pl % get_asFloats ( 'a_nonSchmid' , defaultVal = emptyRealArray )
2020-03-16 16:23:36 +05:30
prm % nonSchmid_pos = lattice_nonSchmidMatrix ( N_sl , a , + 1 )
prm % nonSchmid_neg = lattice_nonSchmidMatrix ( N_sl , a , - 1 )
2019-03-22 15:25:17 +05:30
else
2020-03-16 18:03:14 +05:30
prm % nonSchmid_pos = prm % P_sl
prm % nonSchmid_neg = prm % P_sl
2019-03-22 15:25:17 +05:30
endif
2020-02-14 13:30:14 +05:30
2020-08-15 19:32:10 +05:30
prm % h_sl_sl = lattice_interaction_SlipBySlip ( N_sl , pl % get_asFloats ( 'h_sl_sl' ) , &
phase % get_asString ( 'lattice' ) )
prm % forestProjection = lattice_forestProjection_edge ( N_sl , phase % get_asString ( 'lattice' ) , &
phase % get_asFloat ( 'c/a' , defaultVal = 0.0_pReal ) )
2020-03-16 16:23:36 +05:30
prm % forestProjection = transpose ( prm % forestProjection )
2020-02-14 13:30:14 +05:30
2020-08-15 19:32:10 +05:30
rho_mob_0 = pl % get_asFloats ( 'rho_mob_0' , requiredSize = size ( N_sl ) )
rho_dip_0 = pl % get_asFloats ( 'rho_dip_0' , requiredSize = size ( N_sl ) )
2020-09-22 18:05:05 +05:30
prm % v_0 = pl % get_asFloats ( 'v_0' , requiredSize = size ( N_sl ) )
2020-08-15 19:32:10 +05:30
prm % b_sl = pl % get_asFloats ( 'b_sl' , requiredSize = size ( N_sl ) )
2020-09-22 18:05:05 +05:30
prm % Q_s = pl % get_asFloats ( 'Q_s' , requiredSize = size ( N_sl ) )
2020-03-16 03:37:23 +05:30
2020-08-15 19:32:10 +05:30
prm % i_sl = pl % get_asFloats ( 'i_sl' , requiredSize = size ( N_sl ) )
2020-09-23 04:25:19 +05:30
prm % tau_Peierls = pl % get_asFloats ( 'tau_Peierls' , requiredSize = size ( N_sl ) )
2020-08-15 19:32:10 +05:30
prm % p = pl % get_asFloats ( 'p_sl' , requiredSize = size ( N_sl ) , &
2020-03-16 03:37:23 +05:30
defaultVal = [ ( 1.0_pReal , i = 1 , size ( N_sl ) ) ] )
2020-08-15 19:32:10 +05:30
prm % q = pl % get_asFloats ( 'q_sl' , requiredSize = size ( N_sl ) , &
2020-03-16 03:37:23 +05:30
defaultVal = [ ( 1.0_pReal , i = 1 , size ( N_sl ) ) ] )
2020-09-22 18:05:05 +05:30
prm % h = pl % get_asFloats ( 'h' , requiredSize = size ( N_sl ) )
2020-08-15 19:32:10 +05:30
prm % w = pl % get_asFloats ( 'w' , requiredSize = size ( N_sl ) )
prm % omega = pl % get_asFloats ( 'omega' , requiredSize = size ( N_sl ) )
prm % B = pl % get_asFloats ( 'B' , requiredSize = size ( N_sl ) )
prm % D = pl % get_asFloat ( 'D' )
prm % D_0 = pl % get_asFloat ( 'D_0' )
prm % Q_cl = pl % get_asFloat ( 'Q_cl' )
2020-09-22 18:05:05 +05:30
prm % f_at = pl % get_asFloat ( 'f_at' ) * prm % b_sl ** 3.0_pReal
2020-08-15 19:32:10 +05:30
prm % D_a = pl % get_asFloat ( 'D_a' ) * prm % b_sl
2020-09-14 00:30:34 +05:30
2020-11-29 16:28:39 +05:30
prm % dipoleformation = . not . pl % get_asBool ( 'no_dipole_formation' , defaultVal = . false . )
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
! expand: family => system
2020-03-16 16:23:36 +05:30
rho_mob_0 = math_expand ( rho_mob_0 , N_sl )
rho_dip_0 = math_expand ( rho_dip_0 , N_sl )
2020-03-16 03:37:23 +05:30
prm % q = math_expand ( prm % q , N_sl )
prm % p = math_expand ( prm % p , N_sl )
2020-09-22 18:05:05 +05:30
prm % Q_s = math_expand ( prm % Q_s , N_sl )
2020-03-16 03:37:23 +05:30
prm % b_sl = math_expand ( prm % b_sl , N_sl )
2020-09-22 18:05:05 +05:30
prm % h = math_expand ( prm % h , N_sl )
2020-03-16 03:37:23 +05:30
prm % w = math_expand ( prm % w , N_sl )
prm % omega = math_expand ( prm % omega , N_sl )
2020-09-23 04:25:19 +05:30
prm % tau_Peierls = math_expand ( prm % tau_Peierls , N_sl )
2020-09-22 18:05:05 +05:30
prm % v_0 = math_expand ( prm % v_0 , N_sl )
2020-03-16 03:37:23 +05:30
prm % B = math_expand ( prm % B , N_sl )
prm % i_sl = math_expand ( prm % i_sl , N_sl )
2020-09-22 18:05:05 +05:30
prm % f_at = math_expand ( prm % f_at , N_sl )
2020-03-16 03:37:23 +05:30
prm % D_a = math_expand ( prm % D_a , N_sl )
2020-02-14 13:30:14 +05:30
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'
2020-08-15 19:32:10 +05:30
if ( any ( prm % D_a < = 0.0_pReal ) ) extmsg = trim ( extmsg ) / / ' D_a or 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'
2020-02-14 13:30:14 +05:30
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
2020-09-23 04:25:19 +05:30
allocate ( prm % b_sl , prm % D_a , 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
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2021-03-05 01:46:36 +05:30
Nmembers = count ( material_phaseAt2 == ph )
2019-03-22 15:25:17 +05:30
sizeDotState = size ( [ 'rho_mob ' , 'rho_dip ' , 'gamma_sl' ] ) * prm % sum_N_sl
sizeState = sizeDotState
2020-02-14 13:30:14 +05:30
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
2019-01-06 04:03:18 +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'
2020-02-14 13:30:14 +05:30
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 )
2020-02-14 13:30:14 +05:30
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 ) = 1.0e-2_pReal
2019-03-22 15:25:17 +05:30
! global alias
2021-02-16 20:36:09 +05:30
plasticState ( ph ) % slipRate = > plasticState ( ph ) % dotState ( startIndex : endIndex , : )
2020-02-14 13:30:14 +05:30
2021-03-05 01:46:36 +05:30
allocate ( dst % Lambda_sl ( prm % sum_N_sl , Nmembers ) , source = 0.0_pReal )
allocate ( dst % threshold_stress ( prm % sum_N_sl , Nmembers ) , source = 0.0_pReal )
2020-02-14 13:30:14 +05:30
2021-02-16 20:36:09 +05:30
plasticState ( ph ) % state0 = plasticState ( ph ) % state ! ToDo: this could be done centrally
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
end associate
2020-02-14 13:30:14 +05:30
2020-03-15 00:33:57 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
2020-11-29 16:28:39 +05:30
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
2020-11-29 16:28:39 +05:30
end function plastic_dislotungsten_init
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +05:30
!> @brief Calculate plastic velocity gradient and its tangent.
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2021-01-26 12:03:04 +05:30
pure module subroutine dislotungsten_LpAndItsTangent ( Lp , dLp_dMp , &
2021-01-27 03:11:41 +05:30
Mp , T , ph , me )
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
2020-02-14 13:30:14 +05:30
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 ) :: &
2021-01-27 03:11:41 +05:30
ph , &
2021-01-26 13:09:17 +05:30
me
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
integer :: &
i , k , l , m , n
2021-02-14 11:36:14 +05:30
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
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
associate ( prm = > param ( ph ) )
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
call kinetics ( Mp , T , ph , me , 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 ) &
2020-03-16 18:03:14 +05:30
+ ddot_gamma_dtau_pos ( i ) * prm % P_sl ( k , l , i ) * prm % nonSchmid_pos ( m , n , i ) &
+ ddot_gamma_dtau_neg ( i ) * prm % P_sl ( k , l , i ) * prm % nonSchmid_neg ( m , n , i )
2019-03-22 15:25:17 +05:30
enddo
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
end associate
2018-11-29 03:08:14 +05:30
2021-01-26 12:03:04 +05:30
end subroutine dislotungsten_LpAndItsTangent
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +05:30
!> @brief Calculate the rate of change of microstructure.
2014-11-06 00:41:09 +05:30
!--------------------------------------------------------------------------------------------------
2021-01-27 03:11:41 +05:30
module subroutine dislotungsten_dotState ( Mp , T , ph , me )
2020-02-14 13:30:14 +05:30
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 ) :: &
2021-01-27 03:11:41 +05:30
ph , &
2021-01-26 13:09:17 +05:30
me
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
real ( pReal ) :: &
VacancyDiffusion
2021-02-14 11:36:14 +05:30
real ( pReal ) , dimension ( param ( ph ) % sum_N_sl ) :: &
2019-03-22 15:25:17 +05:30
gdot_pos , gdot_neg , &
tau_pos , &
tau_neg , &
2019-03-22 15:36:08 +05:30
v_cl , &
dot_rho_dip_formation , &
dot_rho_dip_climb , &
dip_distance
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
associate ( prm = > param ( ph ) , stt = > state ( ph ) , &
dot = > dotState ( ph ) , dst = > dependentState ( ph ) )
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
call kinetics ( Mp , T , ph , me , &
2019-03-22 15:25:17 +05:30
gdot_pos , gdot_neg , &
tau_pos_out = tau_pos , tau_neg_out = tau_neg )
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
dot % gamma_sl ( : , me ) = ( gdot_pos + gdot_neg ) ! ToDo: needs to be abs
2019-03-22 15:25:17 +05:30
VacancyDiffusion = prm % D_0 * exp ( - prm % Q_cl / ( kB * T ) )
2020-02-14 13:30:14 +05:30
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
2019-03-22 15:36:08 +05:30
dip_distance = math_clip ( 3.0_pReal * prm % mu * prm % b_sl / ( 1 6.0_pReal * PI * abs ( tau_pos ) ) , &
2020-03-16 04:42:18 +05:30
prm % D_a , & ! lower limit
2021-01-26 13:09:17 +05:30
dst % Lambda_sl ( : , me ) ) ! upper limit
dot_rho_dip_formation = merge ( 2.0_pReal * dip_distance * stt % rho_mob ( : , me ) * abs ( dot % gamma_sl ( : , me ) ) / prm % b_sl , & ! ToDo: ignore region of spontaneous annihilation
2019-03-22 15:36:08 +05:30
0.0_pReal , &
prm % dipoleformation )
2020-09-22 18:05:05 +05:30
v_cl = ( 3.0_pReal * prm % mu * VacancyDiffusion * prm % f_at / ( 2.0_pReal * pi * kB * T ) ) &
2019-03-22 15:36:08 +05:30
* ( 1.0_pReal / ( dip_distance + prm % D_a ) )
2021-01-26 13:09:17 +05:30
dot_rho_dip_climb = ( 4.0_pReal * v_cl * stt % rho_dip ( : , me ) ) / ( dip_distance - prm % D_a ) ! ToDo: Discuss with Franz: Stress dependency?
2019-03-22 15:25:17 +05:30
end where
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
dot % rho_mob ( : , me ) = abs ( dot % gamma_sl ( : , me ) ) / ( prm % b_sl * dst % Lambda_sl ( : , me ) ) & ! multiplication
2019-03-22 15:36:08 +05:30
- dot_rho_dip_formation &
2021-01-26 13:09:17 +05:30
- ( 2.0_pReal * prm % D_a ) / prm % b_sl * stt % rho_mob ( : , me ) * abs ( dot % gamma_sl ( : , me ) ) ! Spontaneous annihilation of 2 single edge dislocations
dot % rho_dip ( : , me ) = dot_rho_dip_formation &
- ( 2.0_pReal * prm % D_a ) / prm % b_sl * stt % rho_dip ( : , me ) * abs ( dot % gamma_sl ( : , me ) ) & ! Spontaneous annihilation of a single edge dislocation with a dipole constituent
2019-03-22 15:36:08 +05:30
- dot_rho_dip_climb
2020-02-14 13:30:14 +05:30
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
2014-11-06 00:41:09 +05:30
2018-12-05 03:00:07 +05:30
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +05:30
!> @brief Calculate derived quantities from state.
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
2021-02-14 11:36:14 +05:30
module subroutine dislotungsten_dependentState ( ph , me )
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
integer , intent ( in ) :: &
2021-02-14 11:36:14 +05:30
ph , &
2021-01-26 13:09:17 +05:30
me
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
real ( pReal ) , dimension ( param ( ph ) % sum_N_sl ) :: &
2019-03-22 15:25:17 +05:30
dislocationSpacing
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
associate ( prm = > param ( ph ) , stt = > state ( ph ) , dst = > dependentState ( ph ) )
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
dislocationSpacing = sqrt ( matmul ( prm % forestProjection , stt % rho_mob ( : , me ) + stt % rho_dip ( : , me ) ) )
dst % threshold_stress ( : , me ) = prm % mu * prm % b_sl &
* sqrt ( matmul ( prm % h_sl_sl , stt % rho_mob ( : , me ) + stt % rho_dip ( : , me ) ) )
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
dst % Lambda_sl ( : , me ) = prm % D / ( 1.0_pReal + prm % D * dislocationSpacing / prm % i_sl )
2020-02-14 13:30:14 +05:30
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
2019-03-10 01:13:31 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +05:30
!> @brief Write results to HDF5 output file.
2019-03-10 01:13:31 +05:30
!--------------------------------------------------------------------------------------------------
2021-02-14 11:36:14 +05:30
module subroutine plastic_dislotungsten_results ( ph , group )
2019-03-10 01:13:31 +05:30
2021-02-14 11:36:14 +05:30
integer , intent ( in ) :: ph
2019-04-04 11:22:36 +05:30
character ( len = * ) , intent ( in ) :: group
2020-02-14 13:30:14 +05:30
2019-03-10 01:13:31 +05:30
integer :: o
2021-02-14 11:36:14 +05:30
associate ( prm = > param ( ph ) , stt = > state ( ph ) , dst = > dependentState ( ph ) )
2020-02-14 13:30:14 +05:30
outputsLoop : do o = 1 , size ( prm % output )
select case ( trim ( prm % output ( o ) ) )
2020-09-14 00:30:34 +05:30
case ( 'rho_mob' )
2020-08-15 19:32:10 +05:30
if ( prm % sum_N_sl > 0 ) call results_writeDataset ( group , stt % rho_mob , trim ( prm % output ( o ) ) , &
2020-02-14 13:30:14 +05:30
'mobile dislocation density' , '1/m²' )
2020-09-14 00:30:34 +05:30
case ( 'rho_dip' )
2020-08-15 19:32:10 +05:30
if ( prm % sum_N_sl > 0 ) call results_writeDataset ( group , stt % rho_dip , trim ( prm % output ( o ) ) , &
2020-02-14 13:30:14 +05:30
'dislocation dipole density' '1/m²' )
2020-09-14 00:30:34 +05:30
case ( 'gamma_sl' )
2020-08-15 19:32:10 +05:30
if ( prm % sum_N_sl > 0 ) call results_writeDataset ( group , stt % gamma_sl , trim ( prm % output ( o ) ) , &
2020-02-14 13:30:14 +05:30
'plastic shear' , '1' )
2020-09-14 00:30:34 +05:30
case ( 'Lambda_sl' )
2020-08-15 19:32:10 +05:30
if ( prm % sum_N_sl > 0 ) call results_writeDataset ( group , dst % Lambda_sl , trim ( prm % output ( o ) ) , &
2020-02-14 13:30:14 +05:30
'mean free path for slip' , 'm' )
2020-09-14 00:30:34 +05:30
case ( 'tau_pass' )
2020-08-15 19:32:10 +05:30
if ( prm % sum_N_sl > 0 ) call results_writeDataset ( group , dst % threshold_stress , trim ( prm % output ( o ) ) , &
2020-02-14 13:30:14 +05:30
'threshold stress for slip' , 'Pa' )
2019-03-10 01:13:31 +05:30
end select
enddo outputsLoop
end associate
2020-11-29 16:28:39 +05:30
end subroutine plastic_dislotungsten_results
2019-03-10 01:13:31 +05:30
2019-01-07 11:54:02 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:30:14 +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
2018-11-27 22:55:06 +05:30
!--------------------------------------------------------------------------------------------------
2021-02-14 11:36:14 +05:30
pure subroutine kinetics ( Mp , T , ph , me , &
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 )
2020-02-14 13:30:14 +05:30
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 ) :: &
2021-02-14 11:36:14 +05:30
ph , &
2021-01-26 13:09:17 +05:30
me
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
real ( pReal ) , intent ( out ) , dimension ( param ( ph ) % sum_N_sl ) :: &
2019-03-22 15:25:17 +05:30
dot_gamma_pos , &
dot_gamma_neg
2021-02-14 11:36:14 +05:30
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
2021-02-14 11:36:14 +05:30
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 , &
t_n , t_k , dtk , dtn , &
needsGoodName ! ToDo: @Karo: any idea?
integer :: j
2020-02-14 13:30:14 +05:30
2021-02-14 11:36:14 +05:30
associate ( prm = > param ( ph ) , stt = > state ( ph ) , dst = > dependentState ( ph ) )
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
do j = 1 , prm % sum_N_sl
2020-03-16 23:13:04 +05:30
tau_pos ( j ) = math_tensordot ( Mp , prm % nonSchmid_pos ( 1 : 3 , 1 : 3 , j ) )
tau_neg ( j ) = math_tensordot ( Mp , prm % nonSchmid_neg ( 1 : 3 , 1 : 3 , j ) )
2019-03-22 15:25:17 +05:30
enddo
2020-02-14 13:30:14 +05:30
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-02-14 13:30:14 +05:30
2020-09-22 18:05:05 +05:30
associate ( BoltzmannRatio = > prm % Q_s / ( kB * T ) , &
2021-01-26 13:09:17 +05:30
dot_gamma_0 = > stt % rho_mob ( : , me ) * prm % b_sl * prm % v_0 , &
effectiveLength = > dst % Lambda_sl ( : , me ) - prm % w )
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
significantPositiveTau : where ( abs ( tau_pos ) - dst % threshold_stress ( : , me ) > tol_math_check )
StressRatio = ( abs ( tau_pos ) - dst % threshold_stress ( : , me ) ) / prm % tau_Peierls
2019-03-22 15:25:17 +05:30
StressRatio_p = StressRatio ** prm % p
StressRatio_pminus1 = StressRatio ** ( prm % p - 1.0_pReal )
needsGoodName = exp ( - BoltzmannRatio * ( 1 - StressRatio_p ) ** prm % q )
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
t_n = prm % b_sl / ( needsGoodName * 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-02-14 13:30:14 +05:30
2020-09-22 18:05:05 +05:30
vel = prm % h / ( t_n + t_k )
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
if ( present ( ddot_gamma_dtau_pos ) ) then
2021-01-26 13:09:17 +05:30
significantPositiveTau2 : where ( abs ( tau_pos ) - dst % threshold_stress ( : , me ) > 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-02-14 13:30:14 +05:30
2020-09-22 18:05:05 +05:30
dvel = - 1.0_pReal * prm % h * ( dtk + dtn ) / ( t_n + t_k ) ** 2.0_pReal
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
2021-01-26 13:09:17 +05:30
significantNegativeTau : where ( abs ( tau_neg ) - dst % threshold_stress ( : , me ) > tol_math_check )
StressRatio = ( abs ( tau_neg ) - dst % threshold_stress ( : , me ) ) / prm % tau_Peierls
2019-03-22 15:25:17 +05:30
StressRatio_p = StressRatio ** prm % p
StressRatio_pminus1 = StressRatio ** ( prm % p - 1.0_pReal )
needsGoodName = exp ( - BoltzmannRatio * ( 1 - StressRatio_p ) ** prm % q )
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
t_n = prm % b_sl / ( needsGoodName * 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-02-14 13:30:14 +05:30
2020-09-22 18:05:05 +05:30
vel = prm % h / ( t_n + t_k )
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
if ( present ( ddot_gamma_dtau_neg ) ) then
2021-01-26 13:09:17 +05:30
significantNegativeTau2 : where ( abs ( tau_neg ) - dst % threshold_stress ( : , me ) > 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-02-14 13:30:14 +05:30
2020-09-22 18:05:05 +05:30
dvel = - 1.0_pReal * prm % h * ( dtk + dtn ) / ( t_n + t_k ) ** 2.0_pReal
2020-02-14 13:30:14 +05:30
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
2020-02-14 13:30:14 +05:30
2019-03-22 15:25:17 +05:30
end associate
end associate
2018-12-05 02:03:32 +05:30
2018-11-27 22:55:06 +05:30
end subroutine kinetics
2014-11-06 00:41:09 +05:30
2021-01-26 05:41:32 +05:30
end submodule dislotungsten