DAMASK_EICMD/src/phase_mechanical_plastic_di...

1069 lines
51 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Su Leen Wong, Max-Planck-Institut für Eisenforschung GmbH
!> @author Nan Jia, Max-Planck-Institut für Eisenforschung GmbH
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine incoprorating dislocation and twinning physics
!> @details to be done
!--------------------------------------------------------------------------------------------------
submodule(phase:plastic) dislotwin
2020-02-14 13:53:09 +05:30
type :: tParameters
real(pReal) :: &
Q_cl = 1.0_pReal, & !< activation energy for dislocation climb
2020-03-16 19:28:42 +05:30
omega = 1.0_pReal, & !< frequency factor for dislocation climb
D = 1.0_pReal, & !< grain size
p_sb = 1.0_pReal, & !< p-exponent in shear band velocity
q_sb = 1.0_pReal, & !< q-exponent in shear band velocity
2020-08-15 19:32:10 +05:30
i_tw = 1.0_pReal, & !< adjustment parameter to calculate MFP for twinning
2020-03-16 19:28:42 +05:30
L_tw = 1.0_pReal, & !< Length of twin nuclei in Burgers vectors
L_tr = 1.0_pReal, & !< Length of trans nuclei in Burgers vectors
x_c_tw = 1.0_pReal, & !< critical distance for formation of twin nucleus
x_c_tr = 1.0_pReal, & !< critical distance for formation of trans nucleus
2020-03-16 19:28:42 +05:30
V_cs = 1.0_pReal, & !< cross slip volume
xi_sb = 1.0_pReal, & !< value for shearband resistance
v_sb = 1.0_pReal, & !< value for shearband velocity_0
2020-03-16 19:28:42 +05:30
E_sb = 1.0_pReal, & !< activation energy for shear bands
2021-12-31 02:39:53 +05:30
delta_G = 1.0_pReal, & !< free energy difference between austensite and martensite
2020-08-15 19:32:10 +05:30
i_tr = 1.0_pReal, & !< adjustment parameter to calculate MFP for transformation
2021-12-31 02:39:53 +05:30
h = 1.0_pReal, & !< stack height of hex nucleus
2021-11-25 10:51:56 +05:30
T_ref = T_ROOM, &
a_cI = 1.0_pReal, &
a_cF = 1.0_pReal
real(pReal), dimension(2) :: &
Gamma_sf = 0.0_pReal
2020-02-14 13:53:09 +05:30
real(pReal), allocatable, dimension(:) :: &
2020-09-23 05:36:03 +05:30
b_sl, & !< absolute length of Burgers vector [m] for each slip system
b_tw, & !< absolute length of Burgers vector [m] for each twin system
b_tr, & !< absolute length of Burgers vector [m] for each transformation system
2021-07-24 01:06:26 +05:30
Q_sl,& !< activation energy for glide [J] for each slip system
v_0, & !< dislocation velocity prefactor [m/s] for each slip system
dot_N_0_tw, & !< twin nucleation rate [1/m³s] for each twin system
t_tw, & !< twin thickness [m] for each twin system
i_sl, & !< Adj. parameter for distance between 2 forest dislocations for each slip system
t_tr, & !< martensite lamellar thickness [m] for each trans system
p, & !< p-exponent in glide velocity
q, & !< q-exponent in glide velocity
2021-12-31 02:39:53 +05:30
r, & !< exponent in twin nucleation rate
s, & !< exponent in trans nucleation rate
2021-02-24 16:36:46 +05:30
tau_0, & !< strength due to elements in solid solution
gamma_char, & !< characteristic shear for twins
2021-07-24 01:06:26 +05:30
B, & !< drag coefficient
d_caron !< distance of spontaneous annhihilation
2020-02-14 13:53:09 +05:30
real(pReal), allocatable, dimension(:,:) :: &
2020-08-15 19:32:10 +05:30
h_sl_sl, & !< components of slip-slip interaction matrix
h_sl_tw, & !< components of slip-twin interaction matrix
h_tw_tw, & !< components of twin-twin interaction matrix
h_sl_tr, & !< components of slip-trans interaction matrix
h_tr_tr, & !< components of trans-trans interaction matrix
n0_sl, & !< slip system normal
forestProjection
2020-02-14 13:53:09 +05:30
real(pReal), allocatable, dimension(:,:,:) :: &
P_sl, &
P_tw, &
P_tr
2020-02-14 13:53:09 +05:30
integer :: &
sum_N_sl, & !< total number of active slip system
sum_N_tw, & !< total number of active twin system
2020-02-14 13:53:09 +05:30
sum_N_tr !< total number of active transformation system
integer, allocatable, dimension(:) :: &
N_tw, &
N_tr
2020-02-14 13:53:09 +05:30
integer, allocatable, dimension(:,:) :: &
fcc_twinNucleationSlipPair ! ToDo: Better name? Is also use for trans
character(len=:), allocatable :: &
lattice_tr
2020-02-14 13:53:09 +05:30
character(len=pStringLen), allocatable, dimension(:) :: &
output
logical :: &
ExtendedDislocations, & !< consider split into partials for climb calculation
fccTwinTransNucleation, & !< twinning and transformation models are for fcc
omitDipoles !< flag controlling consideration of dipole formation
character(len=:), allocatable, dimension(:) :: &
systems_sl, &
systems_tw
end type !< container type for internal constitutive parameters
2020-02-14 13:53:09 +05:30
type :: tDislotwinState
real(pReal), dimension(:,:), pointer :: &
rho_mob, &
rho_dip, &
gamma_sl, &
f_tw, &
f_tr
end type tDislotwinState
2020-02-14 13:53:09 +05:30
2021-07-27 02:26:40 +05:30
type :: tDislotwinDependentState
real(pReal), dimension(:,:), allocatable :: &
Lambda_sl, & !< mean free path between 2 obstacles seen by a moving dislocation
Lambda_tw, & !< mean free path between 2 obstacles seen by a growing twin
Lambda_tr, & !< mean free path between 2 obstacles seen by a growing martensite
2020-08-15 19:32:10 +05:30
tau_pass, & !< threshold stress for slip
tau_hat_tw, & !< threshold stress for twinning
tau_hat_tr, & !< threshold stress for transformation
V_tw, & !< volume of a new twin
2021-12-31 03:06:06 +05:30
V_tr !< volume of a new martensite disc
2021-07-27 02:26:40 +05:30
end type tDislotwinDependentState
2020-02-14 13:53:09 +05:30
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:) :: param
type(tDislotwinState), allocatable, dimension(:) :: &
dotState, &
state
2021-07-27 02:26:40 +05:30
type(tDislotwinDependentState), allocatable, dimension(:) :: dependentState
contains
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Perform module initialization.
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2020-08-15 19:32:10 +05:30
module function plastic_dislotwin_init() result(myPlasticity)
2014-06-11 17:41:14 +05:30
2020-08-15 19:32:10 +05:30
logical, dimension(:), allocatable :: myPlasticity
integer :: &
2021-02-16 20:36:09 +05:30
ph, i, &
2021-03-05 01:46:36 +05:30
Nmembers, &
sizeState, sizeDotState, &
startIndex, endIndex
2020-03-16 19:28:42 +05:30
integer, dimension(:), allocatable :: &
N_sl
2020-03-16 19:28:42 +05:30
real(pReal), allocatable, dimension(:) :: &
rho_mob_0, & !< initial unipolar dislocation density per slip system
rho_dip_0 !< initial dipole dislocation density per slip system
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:53:09 +05:30
2020-08-15 19:32:10 +05:30
myPlasticity = plastic_active('dislotwin')
2021-11-22 02:19:04 +05:30
if (count(myPlasticity) == 0) return
print'(/,1x,a)', '<<<+- phase:mechanical:plastic:dislotwin init -+>>>'
print'(/,a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
2021-02-13 17:37:35 +05:30
print'(/,1x,a)', 'A. Ma and F. Roters, Acta Materialia 52(12):36033612, 2004'
print'( 1x,a)', 'https://doi.org/10.1016/j.actamat.2004.04.012'
print'(/,1x,a)', 'F. Roters et al., Computational Materials Science 39:9195, 2007'
print'( 1x,a)', 'https://doi.org/10.1016/j.commatsci.2006.04.014'
print'(/,1x,a)', 'S.L. Wong et al., Acta Materialia 118:140151, 2016'
print'( 1x,a)', 'https://doi.org/10.1016/j.actamat.2016.07.032'
2020-02-14 13:53:09 +05:30
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
2021-11-22 02:19:04 +05:30
if (.not. myPlasticity(ph)) cycle
2021-02-16 20:36:09 +05:30
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
2020-03-15 14:09:35 +05:30
2020-02-14 13:53:09 +05:30
2018-10-18 02:43:47 +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))
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_sl = lattice_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph))
prm%h_sl_sl = lattice_interaction_SlipBySlip(N_sl,pl%get_as1dFloat('h_sl-sl'),phase_lattice(ph))
prm%forestProjection = lattice_forestProjection_edge(N_sl,phase_lattice(ph),phase_cOverA(ph))
prm%forestProjection = transpose(prm%forestProjection)
2020-02-14 13:53:09 +05:30
2021-07-21 19:16:38 +05:30
prm%n0_sl = lattice_slip_normal(N_sl,phase_lattice(ph),phase_cOverA(ph))
prm%fccTwinTransNucleation = phase_lattice(ph) == 'cF' .and. (N_sl(1) == 12)
2021-07-21 19:16:38 +05:30
if (prm%fccTwinTransNucleation) prm%fcc_twinNucleationSlipPair = lattice_FCC_TWINNUCLEATIONSLIPPAIR
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_sl = pl%get_as1dFloat('Q_sl', requiredSize=size(N_sl))
2021-03-11 22:30:07 +05:30
prm%i_sl = pl%get_as1dFloat('i_sl', requiredSize=size(N_sl))
prm%p = pl%get_as1dFloat('p_sl', requiredSize=size(N_sl))
prm%q = pl%get_as1dFloat('q_sl', requiredSize=size(N_sl))
prm%tau_0 = pl%get_as1dFloat('tau_0', requiredSize=size(N_sl))
prm%B = pl%get_as1dFloat('B', requiredSize=size(N_sl), &
defaultVal=[(0.0_pReal, i=1,size(N_sl))])
2020-02-14 13:53:09 +05:30
2021-07-24 01:06:26 +05:30
prm%Q_cl = pl%get_asFloat('Q_cl')
2020-02-14 13:53:09 +05:30
2021-07-04 02:00:06 +05:30
prm%ExtendedDislocations = pl%get_asBool('extend_dislocations',defaultVal = .false.)
prm%omitDipoles = pl%get_asBool('omit_dipoles',defaultVal = .false.)
2020-03-16 19:28:42 +05:30
! multiplication factor according to crystal structure (nearest neighbors bcc vs fcc/hex)
2020-03-15 14:09:35 +05:30
! details: Argon & Moffat, Acta Metallurgica, Vol. 29, pg 293 to 299, 1981
2020-08-15 19:32:10 +05:30
prm%omega = pl%get_asFloat('omega', defaultVal = 1000.0_pReal) &
* merge(12.0_pReal,8.0_pReal,any(phase_lattice(ph) == ['cF','hP']))
2020-02-14 13:53:09 +05:30
! expand: family => system
2020-03-16 19:28:42 +05:30
rho_mob_0 = math_expand(rho_mob_0, N_sl)
rho_dip_0 = math_expand(rho_dip_0, N_sl)
prm%v_0 = math_expand(prm%v_0, N_sl)
prm%b_sl = math_expand(prm%b_sl, N_sl)
prm%Q_sl = math_expand(prm%Q_sl, N_sl)
prm%i_sl = math_expand(prm%i_sl, N_sl)
prm%p = math_expand(prm%p, N_sl)
prm%q = math_expand(prm%q, N_sl)
2021-02-24 16:36:46 +05:30
prm%tau_0 = math_expand(prm%tau_0, N_sl)
prm%B = math_expand(prm%B, N_sl)
2021-07-24 01:06:26 +05:30
prm%d_caron = pl%get_asFloat('D_a') * prm%b_sl
2020-02-14 13:53:09 +05:30
! sanity checks
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%v_0 < 0.0_pReal)) extmsg = trim(extmsg)//' v_0'
if (any(prm%b_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' b_sl'
2021-07-24 01:06:26 +05:30
if (any(prm%Q_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' Q_sl'
if (any(prm%i_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' i_sl'
if (any(prm%B < 0.0_pReal)) extmsg = trim(extmsg)//' B'
2021-07-24 01:06:26 +05:30
if (any(prm%d_caron < 0.0_pReal)) extmsg = trim(extmsg)//' d_caron(D_a,b_sl)'
if (any(prm%p<=0.0_pReal .or. prm%p>1.0_pReal)) extmsg = trim(extmsg)//' p_sl'
if (any(prm%q< 1.0_pReal .or. prm%q>2.0_pReal)) extmsg = trim(extmsg)//' q_sl'
else slipActive
2020-03-16 19:28:42 +05:30
rho_mob_0 = emptyRealArray; rho_dip_0 = emptyRealArray
2021-07-24 01:06:26 +05:30
allocate(prm%b_sl,prm%Q_sl,prm%v_0,prm%i_sl,prm%p,prm%q,prm%B,source=emptyRealArray)
2020-03-16 19:28:42 +05:30
allocate(prm%forestProjection(0,0),prm%h_sl_sl(0,0))
end if slipActive
2020-02-14 13:53:09 +05:30
2018-10-18 02:43:47 +05:30
!--------------------------------------------------------------------------------------------------
! twin related parameters
prm%N_tw = pl%get_as1dInt('N_tw', defaultVal=emptyIntArray)
prm%sum_N_tw = sum(abs(prm%N_tw))
2020-03-16 19:28:42 +05:30
twinActive: if (prm%sum_N_tw > 0) then
prm%systems_tw = lattice_labels_twin(prm%N_tw,phase_lattice(ph))
prm%P_tw = lattice_SchmidMatrix_twin(prm%N_tw,phase_lattice(ph),phase_cOverA(ph))
prm%h_tw_tw = lattice_interaction_TwinByTwin(prm%N_tw,pl%get_as1dFloat('h_tw-tw'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
2020-02-14 13:53:09 +05:30
prm%b_tw = pl%get_as1dFloat('b_tw', requiredSize=size(prm%N_tw))
prm%t_tw = pl%get_as1dFloat('t_tw', requiredSize=size(prm%N_tw))
prm%r = pl%get_as1dFloat('p_tw', requiredSize=size(prm%N_tw))
2020-02-14 13:53:09 +05:30
prm%x_c_tw = pl%get_asFloat('x_c_tw')
2020-08-15 19:32:10 +05:30
prm%L_tw = pl%get_asFloat('L_tw')
prm%i_tw = pl%get_asFloat('i_tw')
2020-02-14 13:53:09 +05:30
prm%gamma_char= lattice_characteristicShear_Twin(prm%N_tw,phase_lattice(ph),phase_cOverA(ph))
2020-02-14 13:53:09 +05:30
if (.not. prm%fccTwinTransNucleation) then
2021-03-11 22:30:07 +05:30
prm%dot_N_0_tw = pl%get_as1dFloat('dot_N_0_tw')
prm%dot_N_0_tw = math_expand(prm%dot_N_0_tw,prm%N_tw)
endif
2020-02-14 13:53:09 +05:30
! expand: family => system
prm%b_tw = math_expand(prm%b_tw,prm%N_tw)
prm%t_tw = math_expand(prm%t_tw,prm%N_tw)
prm%r = math_expand(prm%r,prm%N_tw)
2020-02-14 13:53:09 +05:30
! sanity checks
if ( prm%x_c_tw < 0.0_pReal) extmsg = trim(extmsg)//' x_c_tw'
if ( prm%L_tw < 0.0_pReal) extmsg = trim(extmsg)//' L_tw'
if ( prm%i_tw < 0.0_pReal) extmsg = trim(extmsg)//' i_tw'
if (any(prm%b_tw < 0.0_pReal)) extmsg = trim(extmsg)//' b_tw'
if (any(prm%t_tw < 0.0_pReal)) extmsg = trim(extmsg)//' t_tw'
2020-08-15 19:32:10 +05:30
if (any(prm%r < 0.0_pReal)) extmsg = trim(extmsg)//' p_tw'
if (.not. prm%fccTwinTransNucleation) then
if (any(prm%dot_N_0_tw < 0.0_pReal)) extmsg = trim(extmsg)//' dot_N_0_tw'
end if
2020-03-16 19:28:42 +05:30
else twinActive
allocate(prm%gamma_char,prm%b_tw,prm%dot_N_0_tw,prm%t_tw,prm%r,source=emptyRealArray)
allocate(prm%h_tw_tw(0,0))
end if twinActive
2020-02-14 13:53:09 +05:30
2018-10-19 01:50:26 +05:30
!--------------------------------------------------------------------------------------------------
! transformation related parameters
prm%N_tr = pl%get_as1dInt('N_tr', defaultVal=emptyIntArray)
prm%sum_N_tr = sum(abs(prm%N_tr))
2020-03-16 19:28:42 +05:30
transActive: if (prm%sum_N_tr > 0) then
2021-03-11 22:30:07 +05:30
prm%b_tr = pl%get_as1dFloat('b_tr')
prm%b_tr = math_expand(prm%b_tr,prm%N_tr)
2020-02-14 13:53:09 +05:30
2021-12-31 02:39:53 +05:30
prm%h = pl%get_asFloat('h', defaultVal=0.0_pReal) ! ToDo: This is not optional!
prm%i_tr = pl%get_asFloat('i_tr', defaultVal=0.0_pReal) ! ToDo: This is not optional!
prm%delta_G = pl%get_asFloat('delta_G')
2021-12-31 02:39:53 +05:30
prm%x_c_tr = pl%get_asFloat('x_c_tr', defaultVal=0.0_pReal) ! ToDo: This is not optional!
2020-08-15 19:32:10 +05:30
prm%L_tr = pl%get_asFloat('L_tr')
prm%a_cI = pl%get_asFloat('a_cI', defaultVal=0.0_pReal)
prm%a_cF = pl%get_asFloat('a_cF', defaultVal=0.0_pReal)
2020-02-14 13:53:09 +05:30
prm%lattice_tr = pl%get_asString('lattice_tr')
2020-02-14 13:53:09 +05:30
prm%h_tr_tr = lattice_interaction_TransByTrans(prm%N_tr,pl%get_as1dFloat('h_tr-tr'),&
phase_lattice(ph))
2020-02-14 13:53:09 +05:30
prm%P_tr = lattice_SchmidMatrix_trans(prm%N_tr,prm%lattice_tr, &
0.0_pReal, &
prm%a_cI, &
prm%a_cF)
2020-02-14 13:53:09 +05:30
2021-03-11 22:30:07 +05:30
prm%t_tr = pl%get_as1dFloat('t_tr')
prm%t_tr = math_expand(prm%t_tr,prm%N_tr)
2021-03-11 22:30:07 +05:30
prm%s = pl%get_as1dFloat('p_tr',defaultVal=[0.0_pReal])
prm%s = math_expand(prm%s,prm%N_tr)
! sanity checks
if ( prm%x_c_tr < 0.0_pReal) extmsg = trim(extmsg)//' x_c_tr'
if ( prm%L_tr < 0.0_pReal) extmsg = trim(extmsg)//' L_tr'
if ( prm%i_tr < 0.0_pReal) extmsg = trim(extmsg)//' i_tr'
if (any(prm%t_tr < 0.0_pReal)) extmsg = trim(extmsg)//' t_tr'
2020-08-15 19:32:10 +05:30
if (any(prm%s < 0.0_pReal)) extmsg = trim(extmsg)//' p_tr'
2020-03-16 19:28:42 +05:30
else transActive
2021-12-31 02:39:53 +05:30
allocate(prm%s,prm%b_tr,prm%t_tr,source=emptyRealArray)
allocate(prm%h_tr_tr(0,0))
end if transActive
2020-02-14 13:53:09 +05:30
2019-01-29 11:11:27 +05:30
!--------------------------------------------------------------------------------------------------
! shearband related parameters
prm%v_sb = pl%get_asFloat('v_sb',defaultVal=0.0_pReal)
if (prm%v_sb > 0.0_pReal) then
prm%xi_sb = pl%get_asFloat('xi_sb')
2020-08-15 19:32:10 +05:30
prm%E_sb = pl%get_asFloat('Q_sb')
prm%p_sb = pl%get_asFloat('p_sb')
prm%q_sb = pl%get_asFloat('q_sb')
2020-02-14 13:53:09 +05:30
! sanity checks
if (prm%xi_sb < 0.0_pReal) extmsg = trim(extmsg)//' xi_sb'
2020-08-15 19:32:10 +05:30
if (prm%E_sb < 0.0_pReal) extmsg = trim(extmsg)//' Q_sb'
if (prm%p_sb <= 0.0_pReal) extmsg = trim(extmsg)//' p_sb'
if (prm%q_sb <= 0.0_pReal) extmsg = trim(extmsg)//' q_sb'
end if
2020-02-14 13:53:09 +05:30
2020-03-16 19:28:42 +05:30
!--------------------------------------------------------------------------------------------------
! parameters required for several mechanisms and their interactions
2021-11-22 02:19:04 +05:30
if (prm%sum_N_sl + prm%sum_N_tw + prm%sum_N_tw > 0) &
2020-08-15 19:32:10 +05:30
prm%D = pl%get_asFloat('D')
2020-03-16 19:28:42 +05:30
2021-07-04 02:00:06 +05:30
if (prm%sum_N_tw + prm%sum_N_tr > 0) &
prm%V_cs = pl%get_asFloat('V_cs')
if (prm%sum_N_tw + prm%sum_N_tr > 0 .or. prm%ExtendedDislocations) then
2021-07-04 02:20:14 +05:30
prm%T_ref = pl%get_asFloat('T_ref')
prm%Gamma_sf(1) = pl%get_asFloat('Gamma_sf')
prm%Gamma_sf(2) = pl%get_asFloat('Gamma_sf,T',defaultVal=0.0_pReal)
end if
2020-03-16 19:28:42 +05:30
slipAndTwinActive: if (prm%sum_N_sl * prm%sum_N_tw > 0) then
prm%h_sl_tw = lattice_interaction_SlipByTwin(N_sl,prm%N_tw,pl%get_as1dFloat('h_sl-tw'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
if (prm%fccTwinTransNucleation .and. size(prm%N_tw) /= 1) extmsg = trim(extmsg)//' N_tw: nucleation'
2020-03-16 19:28:42 +05:30
endif slipAndTwinActive
slipAndTransActive: if (prm%sum_N_sl * prm%sum_N_tr > 0) then
prm%h_sl_tr = lattice_interaction_SlipByTrans(N_sl,prm%N_tr,pl%get_as1dFloat('h_sl-tr'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
if (prm%fccTwinTransNucleation .and. size(prm%N_tr) /= 1) extmsg = trim(extmsg)//' N_tr: nucleation'
2020-03-16 19:28:42 +05:30
endif slipAndTransActive
2020-02-14 13:53:09 +05:30
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2021-04-06 15:08:44 +05:30
Nmembers = count(material_phaseID == ph)
sizeDotState = size(['rho_mob ','rho_dip ','gamma_sl']) * prm%sum_N_sl &
+ size(['f_tw']) * prm%sum_N_tw &
+ size(['f_tr']) * prm%sum_N_tr
sizeState = sizeDotState
2020-02-14 13:53:09 +05:30
2021-03-05 01:46:36 +05:30
call phase_allocateState(plasticState(ph),Nmembers,sizeState,sizeDotState,0)
2020-02-14 13:53:09 +05:30
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
! locally defined state aliases and initialization of state0 and atol
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:53:09 +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:53:09 +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)
2021-11-22 02:19:04 +05:30
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_gamma'
2020-02-14 13:53:09 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_tw
2021-02-16 20:36:09 +05:30
stt%f_tw=>plasticState(ph)%state(startIndex:endIndex,:)
dot%f_tw=>plasticState(ph)%dotState(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_f_tw',defaultVal=1.0e-6_pReal)
2021-02-27 14:34:08 +05:30
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_f_tw'
2020-02-14 13:53:09 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_tr
2021-02-16 20:36:09 +05:30
stt%f_tr=>plasticState(ph)%state(startIndex:endIndex,:)
dot%f_tr=>plasticState(ph)%dotState(startIndex:endIndex,:)
2021-02-26 11:39:04 +05:30
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_f_tr',defaultVal=1.0e-6_pReal)
2021-02-27 14:34:08 +05:30
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_f_tr'
2020-02-14 13:53:09 +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%tau_pass (prm%sum_N_sl,Nmembers),source=0.0_pReal)
2020-02-14 13:53:09 +05:30
2021-03-05 01:46:36 +05:30
allocate(dst%Lambda_tw (prm%sum_N_tw,Nmembers),source=0.0_pReal)
allocate(dst%tau_hat_tw (prm%sum_N_tw,Nmembers),source=0.0_pReal)
allocate(dst%V_tw (prm%sum_N_tw,Nmembers),source=0.0_pReal)
2020-02-14 13:53:09 +05:30
2021-03-05 01:46:36 +05:30
allocate(dst%Lambda_tr (prm%sum_N_tr,Nmembers),source=0.0_pReal)
allocate(dst%tau_hat_tr (prm%sum_N_tr,Nmembers),source=0.0_pReal)
allocate(dst%V_tr (prm%sum_N_tr,Nmembers),source=0.0_pReal)
2020-02-14 13:53:09 +05:30
end associate
2020-02-14 13:53:09 +05:30
2020-03-15 14:09:35 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
2020-08-15 19:32:10 +05:30
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(dislotwin)')
2020-03-15 14:09:35 +05:30
end do
2019-01-27 16:07:50 +05:30
2020-08-15 19:32:10 +05:30
end function plastic_dislotwin_init
2014-07-22 13:13:03 +05:30
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Return the homogenized elasticity matrix.
!--------------------------------------------------------------------------------------------------
2021-04-13 00:51:15 +05:30
module function plastic_dislotwin_homogenizedC(ph,en) result(homogenizedC)
2020-02-14 13:53:09 +05:30
integer, intent(in) :: &
2021-04-13 00:51:15 +05:30
ph, en
real(pReal), dimension(6,6) :: &
2021-11-22 02:19:04 +05:30
homogenizedC, &
C
real(pReal), dimension(:,:,:), allocatable :: &
C66_tw, &
C66_tr
integer :: i
2021-12-20 03:23:48 +05:30
real(pReal) :: f_matrix
2020-02-14 13:53:09 +05:30
C = elastic_C66(ph,en)
associate(prm => param(ph), stt => state(ph))
2021-12-20 03:23:48 +05:30
f_matrix = 1.0_pReal &
- sum(stt%f_tw(1:prm%sum_N_tw,en)) &
- sum(stt%f_tr(1:prm%sum_N_tr,en))
2020-02-14 13:53:09 +05:30
2021-12-20 03:23:48 +05:30
homogenizedC = f_matrix * C
twinActive: if (prm%sum_N_tw > 0) then
C66_tw = lattice_C66_twin(prm%N_tw,C,phase_lattice(ph),phase_cOverA(ph))
do i=1,prm%sum_N_tw
homogenizedC = homogenizedC &
+ stt%f_tw(i,en)*C66_tw(1:6,1:6,i)
end do
2021-11-22 02:19:04 +05:30
end if twinActive
transActive: if (prm%sum_N_tr > 0) then
C66_tr = lattice_C66_trans(prm%N_tr,C,prm%lattice_tr,0.0_pReal,prm%a_cI,prm%a_cF)
do i=1,prm%sum_N_tr
homogenizedC = homogenizedC &
+ stt%f_tr(i,en)*C66_tr(1:6,1:6,i)
end do
end if transActive
2020-02-14 13:53:09 +05:30
end associate
2020-02-14 13:53:09 +05:30
2019-01-27 16:45:11 +05:30
end function plastic_dislotwin_homogenizedC
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
module subroutine dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
2020-02-14 13:53:09 +05:30
real(pReal), dimension(3,3), intent(out) :: Lp
real(pReal), dimension(3,3,3,3), intent(out) :: dLp_dMp
real(pReal), dimension(3,3), intent(in) :: Mp
2021-04-25 11:36:52 +05:30
integer, intent(in) :: ph,en
2020-02-14 13:53:09 +05:30
integer :: i,k,l,m,n
real(pReal) :: &
2021-12-20 03:23:48 +05:30
f_matrix,StressRatio_p,&
E_kB_T, &
ddot_gamma_dtau, &
tau, &
T
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 18:22:09 +05:30
dot_gamma_sl,ddot_gamma_dtau_sl
real(pReal), dimension(param(ph)%sum_N_tw) :: &
2021-02-27 14:34:08 +05:30
dot_gamma_tw,ddot_gamma_dtau_tw
real(pReal), dimension(param(ph)%sum_N_tr) :: &
2021-02-27 14:34:08 +05:30
dot_gamma_tr,ddot_gamma_dtau_tr
real(pReal):: dot_gamma_sb
real(pReal), dimension(3,3) :: eigVectors, P_sb
real(pReal), dimension(3) :: eigValues
real(pReal), dimension(3,6), parameter :: &
sb_sComposition = &
reshape(real([&
1, 0, 1, &
1, 0,-1, &
1, 1, 0, &
1,-1, 0, &
0, 1, 1, &
0, 1,-1 &
],pReal),[ 3,6]), &
sb_mComposition = &
reshape(real([&
1, 0,-1, &
1, 0,+1, &
1,-1, 0, &
1, 1, 0, &
0, 1,-1, &
0, 1, 1 &
],pReal),[ 3,6])
2020-02-14 13:53:09 +05:30
T = thermal_T(ph,en)
Lp = 0.0_pReal
2020-02-14 13:53:09 +05:30
dLp_dMp = 0.0_pReal
associate(prm => param(ph), stt => state(ph))
2020-02-14 13:53:09 +05:30
2021-12-20 03:23:48 +05:30
f_matrix = 1.0_pReal &
- sum(stt%f_tw(1:prm%sum_N_tw,en)) &
- sum(stt%f_tr(1:prm%sum_N_tr,en))
call kinetics_sl(Mp,T,ph,en,dot_gamma_sl,ddot_gamma_dtau_sl)
slipContribution: do i = 1, prm%sum_N_sl
Lp = Lp + dot_gamma_sl(i)*prm%P_sl(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
+ ddot_gamma_dtau_sl(i) * prm%P_sl(k,l,i) * prm%P_sl(m,n,i)
end do slipContribution
call kinetics_tw(Mp,T,dot_gamma_sl,ph,en,dot_gamma_tw,ddot_gamma_dtau_tw)
twinContibution: do i = 1, prm%sum_N_tw
Lp = Lp + dot_gamma_tw(i)*prm%P_tw(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
+ ddot_gamma_dtau_tw(i)* prm%P_tw(k,l,i)*prm%P_tw(m,n,i)
end do twinContibution
call kinetics_tr(Mp,T,dot_gamma_sl,ph,en,dot_gamma_tr,ddot_gamma_dtau_tr)
transContibution: do i = 1, prm%sum_N_tr
Lp = Lp + dot_gamma_tr(i)*prm%P_tr(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
+ ddot_gamma_dtau_tr(i)* prm%P_tr(k,l,i)*prm%P_tr(m,n,i)
end do transContibution
2021-12-20 03:23:48 +05:30
Lp = Lp * f_matrix
dLp_dMp = dLp_dMp * f_matrix
shearBandingContribution: if (dNeq0(prm%v_sb)) then
E_kB_T = prm%E_sb/(K_B*T)
call math_eigh33(eigValues,eigVectors,Mp) ! is Mp symmetric by design?
do i = 1,6
P_sb = 0.5_pReal * math_outer(matmul(eigVectors,sb_sComposition(1:3,i)),&
matmul(eigVectors,sb_mComposition(1:3,i)))
tau = math_tensordot(Mp,P_sb)
significantShearBandStress: if (abs(tau) > tol_math_check) then
StressRatio_p = (abs(tau)/prm%xi_sb)**prm%p_sb
dot_gamma_sb = sign(prm%v_sb*exp(-E_kB_T*(1-StressRatio_p)**prm%q_sb), tau)
ddot_gamma_dtau = abs(dot_gamma_sb)*E_kB_T*prm%p_sb*prm%q_sb/prm%xi_sb &
* (abs(tau)/prm%xi_sb)**(prm%p_sb-1.0_pReal) &
* (1.0_pReal-StressRatio_p)**(prm%q_sb-1.0_pReal)
Lp = Lp + dot_gamma_sb * P_sb
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 * P_sb(k,l) * P_sb(m,n)
end if significantShearBandStress
end do
end if shearBandingContribution
end associate
2020-02-14 13:53:09 +05:30
2021-01-26 12:03:04 +05:30
end subroutine dislotwin_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate the rate of change of microstructure.
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine dislotwin_dotState(Mp,T,ph,en)
real(pReal), dimension(3,3), intent(in):: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
T !< temperature at integration point
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
integer :: i
real(pReal) :: &
2021-12-20 03:23:48 +05:30
f_matrix, &
2021-07-24 01:06:26 +05:30
d_hat, &
2020-02-14 13:53:09 +05:30
v_cl, & !< climb velocity
tau, &
2020-02-14 13:53:09 +05:30
sigma_cl, & !< climb stress
2020-09-23 05:36:03 +05:30
b_d !< ratio of Burgers vector to stacking fault width
real(pReal), dimension(param(ph)%sum_N_sl) :: &
dot_rho_dip_formation, &
dot_rho_dip_climb, &
dot_gamma_sl
real(pReal), dimension(param(ph)%sum_N_tw) :: &
2021-02-27 14:34:08 +05:30
dot_gamma_tw
real(pReal), dimension(param(ph)%sum_N_tr) :: &
dot_gamma_tr
real(pReal) :: &
mu, &
nu, &
Gamma
2019-03-22 11:18:38 +05:30
2021-07-24 01:06:26 +05:30
associate(prm => param(ph), stt => state(ph), dot => dotState(ph), dst => dependentState(ph))
2018-09-01 14:15:34 +05:30
mu = elastic_mu(ph,en)
nu = elastic_nu(ph,en)
Gamma = prm%Gamma_sf(1) + prm%Gamma_sf(2) * (T-prm%T_ref)
2021-12-20 03:23:48 +05:30
f_matrix = 1.0_pReal &
- sum(stt%f_tw(1:prm%sum_N_tw,en)) &
- sum(stt%f_tr(1:prm%sum_N_tr,en))
call kinetics_sl(Mp,T,ph,en,dot_gamma_sl)
dot%gamma_sl(:,en) = abs(dot_gamma_sl)
slipState: do i = 1, prm%sum_N_sl
tau = math_tensordot(Mp,prm%P_sl(1:3,1:3,i))
significantSlipStress: if (dEq0(tau) .or. prm%omitDipoles) then
dot_rho_dip_formation(i) = 0.0_pReal
dot_rho_dip_climb(i) = 0.0_pReal
else significantSlipStress
d_hat = 3.0_pReal*mu*prm%b_sl(i)/(16.0_pReal*PI*abs(tau))
d_hat = math_clip(d_hat, right = dst%Lambda_sl(i,en))
d_hat = math_clip(d_hat, left = prm%d_caron(i))
dot_rho_dip_formation(i) = 2.0_pReal*(d_hat-prm%d_caron(i))/prm%b_sl(i) &
* stt%rho_mob(i,en)*abs(dot_gamma_sl(i))
if (dEq(d_hat,prm%d_caron(i))) then
dot_rho_dip_climb(i) = 0.0_pReal
else
! Argon & Moffat, Acta Metallurgica, Vol. 29, pg 293 to 299, 1981
sigma_cl = dot_product(prm%n0_sl(1:3,i),matmul(Mp,prm%n0_sl(1:3,i)))
b_d = merge(24.0_pReal*PI*(1.0_pReal - nu)/(2.0_pReal + nu) * Gamma / (mu*prm%b_sl(i)), &
1.0_pReal, &
prm%ExtendedDislocations)
v_cl = 2.0_pReal*prm%omega*b_d**2*exp(-prm%Q_cl/(K_B*T)) &
* (exp(abs(sigma_cl)*prm%b_sl(i)**3/(K_B*T)) - 1.0_pReal)
dot_rho_dip_climb(i) = 4.0_pReal*v_cl*stt%rho_dip(i,en) &
/ (d_hat-prm%d_caron(i))
end if
end if significantSlipStress
end do slipState
dot%rho_mob(:,en) = abs(dot_gamma_sl)/(prm%b_sl*dst%Lambda_sl(:,en)) &
- dot_rho_dip_formation &
- 2.0_pReal*prm%d_caron/prm%b_sl * stt%rho_mob(:,en)*abs(dot_gamma_sl)
dot%rho_dip(:,en) = dot_rho_dip_formation &
- 2.0_pReal*prm%d_caron/prm%b_sl * stt%rho_dip(:,en)*abs(dot_gamma_sl) &
- dot_rho_dip_climb
call kinetics_tw(Mp,T,dot_gamma_sl,ph,en,dot_gamma_tw)
2021-12-20 03:23:48 +05:30
dot%f_tw(:,en) = f_matrix*dot_gamma_tw/prm%gamma_char
call kinetics_tr(Mp,T,dot_gamma_sl,ph,en,dot_gamma_tr)
2021-12-20 03:23:48 +05:30
dot%f_tr(:,en) = f_matrix*dot_gamma_tr
2018-08-31 19:38:01 +05:30
end associate
2020-02-14 13:53:09 +05:30
2021-01-26 16:47:00 +05:30
end subroutine dislotwin_dotState
2018-09-12 16:55:18 +05:30
2019-01-27 13:05:07 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate derived quantities from state.
2019-01-27 13:05:07 +05:30
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine dislotwin_dependentState(T,ph,en)
2019-01-27 13:05:07 +05:30
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pReal), intent(in) :: &
T
2020-02-14 13:53:09 +05:30
real(pReal) :: &
2021-02-27 14:34:08 +05:30
sumf_tw,Gamma,sumf_tr
real(pReal), dimension(param(ph)%sum_N_sl) :: &
2021-02-27 14:41:43 +05:30
inv_lambda_sl
real(pReal), dimension(param(ph)%sum_N_tw) :: &
inv_lambda_tw_tw, & !< 1/mean free distance between 2 twin stacks from different systems seen by a growing twin
f_over_t_tw
real(pReal), dimension(param(ph)%sum_N_tr) :: &
inv_lambda_tr_tr, & !< 1/mean free distance between 2 martensite stacks from different systems seen by a growing martensite
f_over_t_tr
real(pReal), dimension(:), allocatable :: &
x0
real(pReal) :: &
mu, &
nu
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
2020-02-14 13:53:09 +05:30
mu = elastic_mu(ph,en)
nu = elastic_nu(ph,en)
2021-07-23 03:39:51 +05:30
sumf_tw = sum(stt%f_tw(1:prm%sum_N_tw,en))
sumf_tr = sum(stt%f_tr(1:prm%sum_N_tr,en))
2020-02-14 13:53:09 +05:30
Gamma = prm%Gamma_sf(1) + prm%Gamma_sf(2) * (T-prm%T_ref)
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
!* rescaled volume fraction for topology
f_over_t_tw = stt%f_tw(1:prm%sum_N_tw,en)/prm%t_tw ! this is per system ...
f_over_t_tr = sumf_tr/prm%t_tr ! but this not
2020-03-16 19:28:42 +05:30
! ToDo ...Physically correct, but naming could be adjusted
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
inv_lambda_sl = sqrt(matmul(prm%forestProjection,stt%rho_mob(:,en)+stt%rho_dip(:,en)))/prm%i_sl
if (prm%sum_N_tw > 0 .and. prm%sum_N_sl > 0) &
inv_lambda_sl = inv_lambda_sl + matmul(prm%h_sl_tw,f_over_t_tw)/(1.0_pReal-sumf_tw)
if (prm%sum_N_tr > 0 .and. prm%sum_N_sl > 0) &
inv_lambda_sl = inv_lambda_sl + matmul(prm%h_sl_tr,f_over_t_tr)/(1.0_pReal-sumf_tr)
dst%Lambda_sl(:,en) = prm%D / (1.0_pReal+prm%D*inv_lambda_sl)
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
inv_lambda_tw_tw = matmul(prm%h_tw_tw,f_over_t_tw)/(1.0_pReal-sumf_tw)
dst%Lambda_tw(:,en) = prm%i_tw*prm%D/(1.0_pReal+prm%D*inv_lambda_tw_tw)
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
inv_lambda_tr_tr = matmul(prm%h_tr_tr,f_over_t_tr)/(1.0_pReal-sumf_tr)
dst%Lambda_tr(:,en) = prm%i_tr*prm%D/(1.0_pReal+prm%D*inv_lambda_tr_tr)
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
!* threshold stress for dislocation motion
dst%tau_pass(:,en) = mu*prm%b_sl* sqrt(matmul(prm%h_sl_sl,stt%rho_mob(:,en)+stt%rho_dip(:,en)))
2020-02-14 13:53:09 +05:30
2021-07-23 03:39:51 +05:30
!* threshold stress for growing twin/martensite
dst%tau_hat_tw(:,en) = Gamma/(3.0_pReal*prm%b_tw) &
+ 3.0_pReal*prm%b_tw*mu/(prm%L_tw*prm%b_tw)
dst%tau_hat_tr(:,en) = Gamma/(3.0_pReal*prm%b_tr) &
+ 3.0_pReal*prm%b_tr*mu/(prm%L_tr*prm%b_tr) &
+ prm%h*prm%delta_G/(3.0_pReal*prm%b_tr)
2020-02-14 13:53:09 +05:30
2021-12-31 02:43:48 +05:30
dst%V_tw(:,en) = PI/4.0_pReal*dst%Lambda_tw(:,en)**2*prm%t_tw
dst%V_tr(:,en) = PI/4.0_pReal*dst%Lambda_tr(:,en)**2*prm%t_tr
2020-02-14 13:53:09 +05:30
end associate
2019-01-28 00:14:53 +05:30
2021-01-26 16:47:00 +05:30
end subroutine dislotwin_dependentState
2019-01-27 13:05:07 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_dislotwin_results(ph,group)
integer, intent(in) :: ph
2019-12-03 02:21:25 +05:30
character(len=*), intent(in) :: group
2020-02-14 13:53:09 +05:30
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')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%rho_mob,group,trim(prm%output(ou)), &
'mobile dislocation density','1/m²',prm%systems_sl)
case('rho_dip')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%rho_dip,group,trim(prm%output(ou)), &
'dislocation dipole density','1/m²',prm%systems_sl)
case('gamma_sl')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%gamma_sl,group,trim(prm%output(ou)), &
'plastic shear','1',prm%systems_sl)
case('Lambda_sl')
2021-08-08 01:32:44 +05:30
call results_writeDataset(dst%Lambda_sl,group,trim(prm%output(ou)), &
'mean free path for slip','m',prm%systems_sl)
case('tau_pass')
2021-08-08 01:32:44 +05:30
call results_writeDataset(dst%tau_pass,group,trim(prm%output(ou)), &
'passing stress for slip','Pa',prm%systems_sl)
case('f_tw')
2021-08-08 01:32:44 +05:30
call results_writeDataset(stt%f_tw,group,trim(prm%output(ou)), &
'twinned volume fraction','m³/m³',prm%systems_tw)
case('Lambda_tw')
2021-08-08 01:32:44 +05:30
call results_writeDataset(dst%Lambda_tw,group,trim(prm%output(ou)), &
'mean free path for twinning','m',prm%systems_tw)
case('tau_hat_tw')
2021-08-08 01:32:44 +05:30
call results_writeDataset(dst%tau_hat_tw,group,trim(prm%output(ou)), &
'threshold stress for twinning','Pa',prm%systems_tw)
case('f_tr')
2021-11-22 02:19:04 +05:30
if (prm%sum_N_tr>0) call results_writeDataset(stt%f_tr,group,trim(prm%output(ou)), &
'martensite volume fraction','m³/m³')
end select
end do
end associate
end subroutine plastic_dislotwin_results
2018-09-12 16:55:18 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate shear rates on slip systems, their derivatives with respect to resolved
! stress, and the resolved stress.
2019-01-27 21:29:44 +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-09-12 16:55:18 +05:30
!--------------------------------------------------------------------------------------------------
2021-07-21 18:22:09 +05:30
pure subroutine kinetics_sl(Mp,T,ph,en, &
dot_gamma_sl,ddot_gamma_dtau_sl,tau_sl)
2018-09-12 16:55:18 +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
2020-02-14 13:53:09 +05:30
real(pReal), dimension(param(ph)%sum_N_sl), intent(out) :: &
dot_gamma_sl
real(pReal), dimension(param(ph)%sum_N_sl), optional, intent(out) :: &
2021-07-21 18:22:09 +05:30
ddot_gamma_dtau_sl, &
tau_sl
real(pReal), dimension(param(ph)%sum_N_sl) :: &
ddot_gamma_dtau
2020-02-14 13:53:09 +05:30
real(pReal), dimension(param(ph)%sum_N_sl) :: &
tau, &
stressRatio, &
StressRatio_p, &
Q_kB_T, &
v_wait_inverse, & !< inverse of the effective velocity of a dislocation waiting at obstacles (unsigned)
v_run_inverse, & !< inverse of the velocity of a free moving dislocation (unsigned)
dV_wait_inverse_dTau, &
dV_run_inverse_dTau, &
dV_dTau, &
tau_eff !< effective resolved stress
2020-02-14 13:53:09 +05:30
integer :: i
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
2020-02-14 13:53:09 +05:30
tau = [(math_tensordot(Mp,prm%P_sl(1:3,1:3,i)),i = 1, prm%sum_N_sl)]
2020-02-14 13:53:09 +05:30
tau_eff = abs(tau)-dst%tau_pass(:,en)
2020-02-14 13:53:09 +05:30
significantStress: where(tau_eff > tol_math_check)
stressRatio = tau_eff/prm%tau_0
StressRatio_p = stressRatio** prm%p
Q_kB_T = prm%Q_sl/(K_B*T)
2021-10-31 01:38:55 +05:30
v_wait_inverse = exp(Q_kB_T*(1.0_pReal-StressRatio_p)** prm%q) &
/ prm%v_0
v_run_inverse = prm%B/(tau_eff*prm%b_sl)
2020-02-14 13:53:09 +05:30
dot_gamma_sl = sign(stt%rho_mob(:,en)*prm%b_sl/(v_wait_inverse+v_run_inverse),tau)
2020-02-14 13:53:09 +05:30
dV_wait_inverse_dTau = -1.0_pReal * v_wait_inverse * prm%p * prm%q * Q_kB_T &
* (stressRatio**(prm%p-1.0_pReal)) &
* (1.0_pReal-StressRatio_p)**(prm%q-1.0_pReal) &
/ prm%tau_0
dV_run_inverse_dTau = -1.0_pReal * v_run_inverse/tau_eff
dV_dTau = -1.0_pReal * (dV_wait_inverse_dTau+dV_run_inverse_dTau) &
/ (v_wait_inverse+v_run_inverse)**2
ddot_gamma_dtau = dV_dTau*stt%rho_mob(:,en)*prm%b_sl
else where significantStress
dot_gamma_sl = 0.0_pReal
ddot_gamma_dtau = 0.0_pReal
end where significantStress
2020-02-14 13:53:09 +05:30
end associate
2020-02-14 13:53:09 +05:30
2021-11-22 02:19:04 +05:30
if (present(ddot_gamma_dtau_sl)) ddot_gamma_dtau_sl = ddot_gamma_dtau
if (present(tau_sl)) tau_sl = tau
2020-02-14 13:53:09 +05:30
2021-07-21 18:22:09 +05:30
end subroutine kinetics_sl
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate shear rates on twin 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.
!--------------------------------------------------------------------------------------------------
2021-07-21 18:22:09 +05:30
pure subroutine kinetics_tw(Mp,T,dot_gamma_sl,ph,en,&
dot_gamma_tw,ddot_gamma_dtau_tw)
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), dimension(param(ph)%sum_N_sl), intent(in) :: &
dot_gamma_sl
2020-02-14 13:53:09 +05:30
real(pReal), dimension(param(ph)%sum_N_tw), intent(out) :: &
2021-02-27 14:34:08 +05:30
dot_gamma_tw
real(pReal), dimension(param(ph)%sum_N_tw), optional, intent(out) :: &
2021-02-27 14:34:08 +05:30
ddot_gamma_dtau_tw
2020-02-14 13:53:09 +05:30
real, dimension(param(ph)%sum_N_tw) :: &
tau, &
2021-12-31 02:43:48 +05:30
dot_N_0, &
stressRatio_r, &
ddot_gamma_dtau
2021-12-31 03:06:06 +05:30
real :: &
x0, &
tau_r, &
Gamma, &
mu, nu
2021-12-31 02:43:48 +05:30
integer, dimension(2) :: &
s
integer :: i
2020-02-14 13:53:09 +05:30
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
2020-02-14 13:53:09 +05:30
2021-12-31 03:06:06 +05:30
mu = elastic_mu(ph,en)
nu = elastic_nu(ph,en)
Gamma = prm%Gamma_sf(1) + prm%Gamma_sf(2) * (T-prm%T_ref)
do i = 1, prm%sum_N_tw
tau(i) = math_tensordot(Mp,prm%P_tw(1:3,1:3,i))
isFCC: if (prm%fccTwinTransNucleation) then
2021-12-31 03:06:06 +05:30
x0 = mu*prm%b_tw(i)**2/(Gamma*8.0_pReal*PI)*(2.0_pReal+nu)/(1.0_pReal-nu) ! ToDo: In the paper, this is the Burgers vector for slip
2021-12-31 13:15:55 +05:30
tau_r = mu*prm%b_tw(i)/(2.0_pReal*PI)*(1.0_pReal/(x0+prm%x_c_tw)+cos(PI/3.0_pReal)/x0)
2021-12-31 03:06:06 +05:30
if (tau(i) < tau_r) then ! ToDo: correct?
2021-12-31 02:43:48 +05:30
s=prm%fcc_twinNucleationSlipPair(1:2,i)
dot_N_0(i)=(abs(dot_gamma_sl(s(1)))*(stt%rho_mob(s(2),en)+stt%rho_dip(s(2),en))+&
abs(dot_gamma_sl(s(2)))*(stt%rho_mob(s(1),en)+stt%rho_dip(s(1),en)))/&
(prm%L_tw*prm%b_sl(i))*(1.0_pReal-exp(-prm%V_cs/(K_B*T)*(tau_r-tau(i))))
else
dot_N_0(i)=0.0_pReal
end if
else isFCC
dot_N_0(i)=prm%dot_N_0_tw(i)
end if isFCC
end do
2020-02-14 13:53:09 +05:30
significantStress: where(tau > tol_math_check)
StressRatio_r = (dst%tau_hat_tw(:,en)/tau)**prm%r
2021-12-31 02:43:48 +05:30
dot_gamma_tw = prm%gamma_char * dst%V_tw(:,en) * dot_N_0*exp(-StressRatio_r)
ddot_gamma_dtau = (dot_gamma_tw*prm%r/tau)*StressRatio_r
else where significantStress
dot_gamma_tw = 0.0_pReal
ddot_gamma_dtau = 0.0_pReal
end where significantStress
2020-02-14 13:53:09 +05:30
end associate
2020-02-14 13:53:09 +05:30
2021-11-22 02:19:04 +05:30
if (present(ddot_gamma_dtau_tw)) ddot_gamma_dtau_tw = ddot_gamma_dtau
2018-09-12 16:55:18 +05:30
2021-07-21 18:22:09 +05:30
end subroutine kinetics_tw
2019-01-28 00:14:53 +05:30
2018-09-15 14:13:05 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-14 13:53:09 +05:30
!> @brief Calculate shear rates on transformation 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-09-15 14:13:05 +05:30
!--------------------------------------------------------------------------------------------------
2021-07-21 18:22:09 +05:30
pure subroutine kinetics_tr(Mp,T,dot_gamma_sl,ph,en,&
dot_gamma_tr,ddot_gamma_dtau_tr)
2018-09-15 14:13:05 +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), dimension(param(ph)%sum_N_sl), intent(in) :: &
dot_gamma_sl
2020-02-14 13:53:09 +05:30
real(pReal), dimension(param(ph)%sum_N_tr), intent(out) :: &
dot_gamma_tr
real(pReal), dimension(param(ph)%sum_N_tr), optional, intent(out) :: &
2021-02-27 14:34:08 +05:30
ddot_gamma_dtau_tr
2020-02-14 13:53:09 +05:30
real, dimension(param(ph)%sum_N_tr) :: &
ddot_gamma_dtau
2021-12-31 03:06:06 +05:30
real :: &
2021-12-31 13:15:55 +05:30
stressRatio_s, &
tau, tau_r, &
dot_N_0, &
2021-12-31 03:06:06 +05:30
x0, &
Gamma, &
mu, nu
2021-12-31 02:43:48 +05:30
integer, dimension(2) :: &
s
integer :: i
associate(prm => param(ph), stt => state(ph), dst => dependentState(ph))
2020-02-14 13:53:09 +05:30
2021-12-31 03:06:06 +05:30
mu = elastic_mu(ph,en)
nu = elastic_nu(ph,en)
Gamma = prm%Gamma_sf(1) + prm%Gamma_sf(2) * (T-prm%T_ref)
do i = 1, prm%sum_N_tr
2021-12-31 13:15:55 +05:30
tau = math_tensordot(Mp,prm%P_tr(1:3,1:3,i))
x0 = mu*prm%b_tr(i)**2/(Gamma*8.0_pReal*PI)*(2.0_pReal+nu)/(1.0_pReal-nu) ! ToDo: In the paper, this is the Burgers vector for slip
2021-12-31 13:15:55 +05:30
tau_r = mu*prm%b_tr(i)/(2.0_pReal*PI)*(1.0_pReal/(x0+prm%x_c_tr)+cos(PI/3.0_pReal)/x0)
if (tau > tol_math_check .and. tau < tau_r) then
2021-12-31 02:43:48 +05:30
s=prm%fcc_twinNucleationSlipPair(1:2,i)
2021-12-31 13:15:55 +05:30
dot_N_0=(abs(dot_gamma_sl(s(1)))*(stt%rho_mob(s(2),en)+stt%rho_dip(s(2),en))+&
abs(dot_gamma_sl(s(2)))*(stt%rho_mob(s(1),en)+stt%rho_dip(s(1),en)))/&
(prm%L_tr*prm%b_sl(i))*(1.0_pReal-exp(-prm%V_cs/(K_B*T)*(tau_r-tau)))
StressRatio_s = (dst%tau_hat_tr(i,en)/tau)**prm%s(i)
dot_gamma_tr(i) = dst%V_tr(i,en) * dot_N_0*exp(-StressRatio_s)
ddot_gamma_dtau(i) = (dot_gamma_tr(i)*prm%s(i)/tau)*StressRatio_s
2021-12-31 02:39:53 +05:30
else
2021-12-31 13:15:55 +05:30
dot_gamma_tr(i) = 0.0_pReal
ddot_gamma_dtau(i) = 0.0_pReal
2021-12-31 02:39:53 +05:30
end if
end do
2020-02-14 13:53:09 +05:30
end associate
2020-02-14 13:53:09 +05:30
2021-11-22 02:19:04 +05:30
if (present(ddot_gamma_dtau_tr)) ddot_gamma_dtau_tr = ddot_gamma_dtau
2019-01-28 00:14:53 +05:30
2021-07-21 18:22:09 +05:30
end subroutine kinetics_tr
2021-01-26 05:41:32 +05:30
end submodule dislotwin