DAMASK_EICMD/src/plastic_dislotwin.f90

1322 lines
58 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
!--------------------------------------------------------------------------------------------------
module plastic_dislotwin
use prec, only: &
2019-03-19 02:47:11 +05:30
pReal
2019-01-27 16:07:50 +05:30
implicit none
private
2019-03-19 02:47:11 +05:30
integer, dimension(:,:), allocatable, target, public :: &
2019-01-27 16:07:50 +05:30
plastic_dislotwin_sizePostResult !< size of each post result output
2019-03-19 02:47:11 +05:30
character(len=64), dimension(:,:), allocatable, target, public :: &
2019-01-27 16:07:50 +05:30
plastic_dislotwin_output !< name of each post result output
2019-03-19 02:47:11 +05:30
real(pReal), parameter, private :: &
2019-01-27 16:07:50 +05:30
kB = 1.38e-23_pReal !< Boltzmann constant in J/Kelvin
2018-08-31 19:38:01 +05:30
2019-01-27 16:07:50 +05:30
enum, bind(c)
enumerator :: &
2018-10-18 01:17:50 +05:30
undefined_ID, &
2019-03-19 02:38:41 +05:30
rho_mob_ID, &
rho_dip_ID, &
2019-03-22 11:18:38 +05:30
dot_gamma_sl_ID, &
2019-03-19 02:38:41 +05:30
gamma_sl_ID, &
2019-03-20 14:16:49 +05:30
Lambda_sl_ID, &
2018-10-18 01:17:50 +05:30
resolved_stress_slip_ID, &
threshold_stress_slip_ID, &
edge_dipole_distance_ID, &
2019-03-19 02:38:41 +05:30
f_tw_ID, &
2019-03-20 14:16:49 +05:30
Lambda_tw_ID, &
2018-10-18 01:17:50 +05:30
resolved_stress_twin_ID, &
2019-03-22 11:27:32 +05:30
tau_hat_tw_ID, &
2019-03-20 14:16:49 +05:30
f_tr_ID
end enum
2019-01-27 16:07:50 +05:30
2018-10-18 01:17:50 +05:30
type, private :: tParameters
real(pReal) :: &
mu, &
nu, &
D0, & !< prefactor for self-diffusion coefficient
Qsd, & !< activation energy for dislocation climb
2019-03-21 02:21:52 +05:30
D, & !<grain size
2019-03-21 11:59:38 +05:30
p_sb, & !< p-exponent in shear band velocity
q_sb, & !< q-exponent in shear band velocity
2018-10-18 01:17:50 +05:30
CEdgeDipMinDistance, & !<
2019-03-21 02:21:52 +05:30
i_tw, & !<
2019-03-22 11:27:32 +05:30
tau_0, & !<strength due to elements in solid solution
2019-03-21 11:59:38 +05:30
L_tw, & !< Length of twin nuclei in Burgers vectors
L_tr, & !< Length of trans nuclei in Burgers vectors
2018-10-18 01:17:50 +05:30
xc_twin, & !< critical distance for formation of twin nucleus
xc_trans, & !< critical distance for formation of trans nucleus
2019-03-21 11:59:38 +05:30
V_cs, & !< cross slip volume
2018-10-18 01:17:50 +05:30
sbResistance, & !< value for shearband resistance (might become an internal state variable at some point)
sbVelocity, & !< value for shearband velocity_0
2019-03-21 02:21:52 +05:30
sbQedge, & !< activation energy for shear bands
2018-10-18 01:17:50 +05:30
SFE_0K, & !< stacking fault energy at zero K
dSFE_dT, & !< temperature dependance of stacking fault energy
2019-03-19 10:55:21 +05:30
aTol_rho, & !< absolute tolerance for integration of dislocation density
aTol_f_tw, & !< absolute tolerance for integration of twin volume fraction
aTol_f_tr, & !< absolute tolerance for integration of trans volume fraction
2019-03-21 11:59:38 +05:30
gamma_fcc_hex, & !< Free energy difference between austensite and martensite
2019-03-21 02:21:52 +05:30
i_tr, & !<
2018-10-18 01:17:50 +05:30
transStackHeight !< Stack height of hex nucleus
2019-03-19 02:47:11 +05:30
real(pReal), dimension(:), allocatable :: &
2019-03-19 03:14:54 +05:30
rho_mob_0, & !< initial unipolar dislocation density per slip system
rho_dip_0, & !< initial dipole dislocation density per slip system
b_sl, & !< absolute length of burgers vector [m] for each slip system
b_tw, & !< absolute length of burgers vector [m] for each twin system
2019-03-20 14:16:49 +05:30
b_tr, & !< absolute length of burgers vector [m] for each transformation system
2019-03-21 02:21:52 +05:30
Delta_F,& !< activation energy for glide [J] for each slip system
2018-10-18 01:17:50 +05:30
v0, & !< dislocation velocity prefactor [m/s] for each slip system
Ndot0_twin, & !< twin nucleation rate [1/m³s] for each twin system
Ndot0_trans, & !< trans nucleation rate [1/m³s] for each trans system
twinsize, & !< twin thickness [m] for each twin system
CLambdaSlip, & !< Adj. parameter for distance between 2 forest dislocations for each slip system
atomicVolume, &
lamellarsize, & !< martensite lamellar thickness [m] for each trans system and instance
2018-10-18 01:17:50 +05:30
p, & !< p-exponent in glide velocity
q, & !< q-exponent in glide velocity
r, & !< r-exponent in twin nucleation rate
s, & !< s-exponent in trans nucleation rate
2019-03-21 11:59:38 +05:30
gamma_char, & !< characteristic shear for twins
2018-10-18 01:17:50 +05:30
B !< drag coefficient
2019-03-19 02:47:11 +05:30
real(pReal), dimension(:,:), allocatable :: &
2019-03-19 02:38:41 +05:30
h_sl_sl, & !<
h_sl_tw, & !<
2019-03-21 11:59:38 +05:30
h_tw_tw, & !<
h_sl_tr, & !<
h_tr_tr !<
2019-03-19 02:47:11 +05:30
integer, dimension(:,:), allocatable :: &
fcc_twinNucleationSlipPair ! ToDo: Better name? Is also use for trans
2019-03-19 02:47:11 +05:30
real(pReal), dimension(:,:), allocatable :: &
forestProjection, &
C66
2019-03-19 02:47:11 +05:30
real(pReal), dimension(:,:,:), allocatable :: &
2019-03-21 11:59:38 +05:30
P_tr, &
P_sl, &
P_tw, &
C66_tw, &
C66_tr
2019-03-19 02:47:11 +05:30
integer :: &
2019-03-21 02:21:52 +05:30
sum_N_sl, & !< total number of active slip system
sum_N_tw, & !< total number of active twin system
sum_N_tr !< total number of active transformation system
2019-03-19 02:47:11 +05:30
integer, dimension(:), allocatable :: &
2019-03-19 02:54:45 +05:30
N_sl, & !< number of active slip systems for each family
N_tw, & !< number of active twin systems for each family
2019-03-21 02:21:52 +05:30
N_tr !< number of active transformation systems for each family
2019-03-19 02:47:11 +05:30
integer(kind(undefined_ID)), dimension(:), allocatable :: &
2019-01-27 16:07:50 +05:30
outputID !< ID of each post result output
logical :: &
fccTwinTransNucleation, & !< twinning and transformation models are for fcc
dipoleFormation !< flag indicating consideration of dipole formation
end type !< container type for internal constitutive parameters
type, private :: tDislotwinState
2019-03-19 02:47:11 +05:30
real(pReal), dimension(:,:), pointer :: &
2019-03-21 02:21:52 +05:30
rho_mob, &
rho_dip, &
2019-03-21 11:59:38 +05:30
gamma_sl, &
2019-03-21 02:21:52 +05:30
f_tw, &
2019-03-21 11:59:38 +05:30
f_tr
end type tDislotwinState
type, private :: tDislotwinMicrostructure
2019-03-19 02:47:11 +05:30
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
tau_pass, &
2019-03-22 11:27:32 +05:30
tau_hat_tw, &
tau_hat_tr, &
2019-03-21 02:21:52 +05:30
f_tw, &
f_tr, &
tau_r_tw, & !< stress to bring partials close together (twin)
tau_r_tr !< stress to bring partials close together (trans)
end type tDislotwinMicrostructure
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:), private :: param
2019-03-19 02:47:11 +05:30
type(tDislotwinState), allocatable, dimension(:), private :: &
2019-01-27 16:07:50 +05:30
dotState, &
state
type(tDislotwinMicrostructure), allocatable, dimension(:), private :: microstructure
public :: &
plastic_dislotwin_init, &
plastic_dislotwin_homogenizedC, &
2019-01-27 13:05:07 +05:30
plastic_dislotwin_dependentState, &
plastic_dislotwin_LpAndItsTangent, &
plastic_dislotwin_dotState, &
plastic_dislotwin_postResults, &
plastic_dislotwin_results
2019-01-27 16:07:50 +05:30
private :: &
kinetics_slip, &
kinetics_twin, &
kinetics_trans
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine plastic_dislotwin_init
use prec, only: &
2018-10-18 01:48:33 +05:30
pStringLen, &
dEq0, &
dNeq0, &
dNeq
use debug, only: &
debug_level,&
debug_constitutive,&
debug_levelBasic
use math, only: &
2018-05-17 23:02:41 +05:30
math_expand,&
PI
use IO, only: &
2019-03-09 05:37:26 +05:30
IO_error
use material, only: &
phase_plasticity, &
phase_plasticityInstance, &
phase_Noutput, &
material_allocatePlasticState, &
PLASTICITY_DISLOTWIN_label, &
PLASTICITY_DISLOTWIN_ID, &
2014-06-11 17:41:14 +05:30
material_phase, &
plasticState
use config, only: &
2018-06-27 21:08:52 +05:30
config_phase
use lattice
2014-06-11 17:41:14 +05:30
implicit none
2019-03-19 02:47:11 +05:30
integer :: &
2019-01-27 16:07:50 +05:30
Ninstance, &
p, i, &
NipcMyPhase, outputSize, &
sizeState, sizeDotState, &
startIndex, endIndex
2019-03-19 02:47:11 +05:30
integer, dimension(0), parameter :: emptyIntArray = [integer::]
real(pReal), dimension(0), parameter :: emptyRealArray = [real(pReal)::]
character(len=65536), dimension(0), parameter :: emptyStringArray = [character(len=65536)::]
2018-06-27 21:08:52 +05:30
2018-10-18 01:48:33 +05:30
integer(kind(undefined_ID)) :: &
2019-01-27 16:07:50 +05:30
outputID
2018-10-18 01:48:33 +05:30
character(len=pStringLen) :: &
extmsg = ''
character(len=65536), dimension(:), allocatable :: &
2019-01-27 16:07:50 +05:30
outputs
2018-10-18 01:48:33 +05:30
write(6,'(/,a)') ' <<<+- constitutive_'//PLASTICITY_DISLOTWIN_label//' init -+>>>'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Ma and Roters, Acta Materialia 52(12):36033612, 2004'
2018-07-17 00:44:33 +05:30
write(6,'(a)') ' https://doi.org/10.1016/j.actamat.2004.04.012'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Roters et al., Computational Materials Science 39:9195, 2007'
2018-07-17 00:44:33 +05:30
write(6,'(a)') ' https://doi.org/10.1016/j.commatsci.2006.04.014'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Wong et al., Acta Materialia 118:140151, 2016'
2018-08-31 18:03:42 +05:30
write(6,'(a,/)') ' https://doi.org/10.1016/j.actamat.2016.07.032'
2019-01-27 16:07:50 +05:30
2019-03-09 05:37:26 +05:30
Ninstance = count(phase_plasticity == PLASTICITY_DISLOTWIN_ID)
2019-01-27 16:07:50 +05:30
2019-03-19 02:47:11 +05:30
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
2018-10-18 01:48:33 +05:30
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-03-19 02:47:11 +05:30
allocate(plastic_dislotwin_sizePostResult(maxval(phase_Noutput),Ninstance),source=0)
2018-10-18 01:48:33 +05:30
allocate(plastic_dislotwin_output(maxval(phase_Noutput),Ninstance))
plastic_dislotwin_output = ''
2018-10-18 01:48:33 +05:30
allocate(param(Ninstance))
allocate(state(Ninstance))
allocate(dotState(Ninstance))
allocate(microstructure(Ninstance))
2019-03-19 02:47:11 +05:30
do p = 1, size(phase_plasticity)
2018-06-27 21:08:52 +05:30
if (phase_plasticity(p) /= PLASTICITY_DISLOTWIN_ID) cycle
2018-09-05 19:45:57 +05:30
associate(prm => param(phase_plasticityInstance(p)), &
2018-09-14 15:26:36 +05:30
dot => dotState(phase_plasticityInstance(p)), &
2018-09-05 19:45:57 +05:30
stt => state(phase_plasticityInstance(p)), &
2019-01-28 02:38:36 +05:30
dst => microstructure(phase_plasticityInstance(p)), &
config => config_phase(p))
2019-03-19 10:55:21 +05:30
prm%aTol_rho = config%getFloat('atol_rho', defaultVal=0.0_pReal)
prm%aTol_f_tw = config%getFloat('atol_twinfrac', defaultVal=0.0_pReal)
prm%aTol_f_tr = config%getFloat('atol_transfrac', defaultVal=0.0_pReal)
2019-01-29 11:11:27 +05:30
2018-09-05 19:45:57 +05:30
! This data is read in already in lattice
prm%mu = lattice_mu(p)
prm%nu = lattice_nu(p)
prm%C66 = lattice_C66(1:6,1:6,p)
2018-10-18 02:43:47 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
2019-03-19 02:54:45 +05:30
prm%N_sl = config%getInts('nslip',defaultVal=emptyIntArray)
2019-03-20 14:16:49 +05:30
prm%sum_N_sl = sum(prm%N_sl)
slipActive: if (prm%sum_N_sl > 0) then
2019-03-21 11:59:38 +05:30
prm%P_sl = lattice_SchmidMatrix_slip(prm%N_sl,config%getString('lattice_structure'),&
config%getFloat('c/a',defaultVal=0.0_pReal))
2019-03-19 02:54:45 +05:30
prm%h_sl_sl = lattice_interaction_SlipBySlip(prm%N_sl, &
2019-03-19 02:38:41 +05:30
config%getFloats('interaction_slipslip'), &
config%getString('lattice_structure'))
2019-03-19 02:54:45 +05:30
prm%forestProjection = lattice_forestProjection (prm%N_sl,config%getString('lattice_structure'),&
2019-01-27 16:07:50 +05:30
config%getFloat('c/a',defaultVal=0.0_pReal))
prm%fccTwinTransNucleation = merge(.true., .false., lattice_structure(p) == LATTICE_FCC_ID) &
2019-03-19 02:54:45 +05:30
.and. (prm%N_sl(1) == 12)
2019-01-27 16:07:50 +05:30
if(prm%fccTwinTransNucleation) &
prm%fcc_twinNucleationSlipPair = lattice_fcc_twinNucleationSlipPair
2018-10-18 02:43:47 +05:30
2019-03-19 03:14:54 +05:30
prm%rho_mob_0 = config%getFloats('rhoedge0', requiredSize=size(prm%N_sl))
prm%rho_dip_0 = config%getFloats('rhoedgedip0',requiredSize=size(prm%N_sl))
2019-03-19 02:54:45 +05:30
prm%v0 = config%getFloats('v0', requiredSize=size(prm%N_sl))
2019-03-19 03:14:54 +05:30
prm%b_sl = config%getFloats('slipburgers',requiredSize=size(prm%N_sl))
2019-03-21 02:21:52 +05:30
prm%Delta_F = config%getFloats('qedge', requiredSize=size(prm%N_sl))
2019-03-19 02:54:45 +05:30
prm%CLambdaSlip = config%getFloats('clambdaslip',requiredSize=size(prm%N_sl))
prm%p = config%getFloats('p_slip', requiredSize=size(prm%N_sl))
prm%q = config%getFloats('q_slip', requiredSize=size(prm%N_sl))
prm%B = config%getFloats('b', requiredSize=size(prm%N_sl), &
defaultVal=[(0.0_pReal, i=1,size(prm%N_sl))])
2018-06-27 21:08:52 +05:30
2019-03-22 11:27:32 +05:30
prm%tau_0 = config%getFloat('solidsolutionstrength')
prm%CEdgeDipMinDistance = config%getFloat('cedgedipmindistance')
2019-01-29 11:11:27 +05:30
prm%D0 = config%getFloat('d0')
prm%Qsd = config%getFloat('qsd')
2019-03-19 03:14:54 +05:30
prm%atomicVolume = config%getFloat('catomicvolume') * prm%b_sl**3.0_pReal
2018-06-27 21:08:52 +05:30
2018-10-18 03:26:57 +05:30
! expand: family => system
2019-03-19 03:14:54 +05:30
prm%rho_mob_0 = math_expand(prm%rho_mob_0, prm%N_sl)
prm%rho_dip_0 = math_expand(prm%rho_dip_0, prm%N_sl)
2019-03-19 02:54:45 +05:30
prm%v0 = math_expand(prm%v0, prm%N_sl)
2019-03-19 03:14:54 +05:30
prm%b_sl = math_expand(prm%b_sl,prm%N_sl)
2019-03-21 02:21:52 +05:30
prm%Delta_F = math_expand(prm%Delta_F, prm%N_sl)
2019-03-19 02:54:45 +05:30
prm%CLambdaSlip = math_expand(prm%CLambdaSlip, prm%N_sl)
prm%p = math_expand(prm%p, prm%N_sl)
prm%q = math_expand(prm%q, prm%N_sl)
prm%B = math_expand(prm%B, prm%N_sl)
prm%atomicVolume = math_expand(prm%atomicVolume,prm%N_sl)
2018-10-12 20:54:46 +05:30
2018-10-18 03:26:57 +05:30
! sanity checks
2019-01-29 11:11:27 +05:30
if ( prm%D0 <= 0.0_pReal) extmsg = trim(extmsg)//' D0'
if ( prm%Qsd <= 0.0_pReal) extmsg = trim(extmsg)//' Qsd'
2019-03-19 03:14:54 +05:30
if (any(prm%rho_mob_0 < 0.0_pReal)) extmsg = trim(extmsg)//' rho_mob_0'
if (any(prm%rho_dip_0 < 0.0_pReal)) extmsg = trim(extmsg)//' rho_dip_0'
2019-01-29 11:11:27 +05:30
if (any(prm%v0 < 0.0_pReal)) extmsg = trim(extmsg)//' v0'
2019-03-19 03:14:54 +05:30
if (any(prm%b_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' b_sl'
2019-03-21 02:21:52 +05:30
if (any(prm%Delta_F <= 0.0_pReal)) extmsg = trim(extmsg)//' Delta_F'
2019-01-29 11:11:27 +05:30
if (any(prm%CLambdaSlip <= 0.0_pReal)) extmsg = trim(extmsg)//' CLambdaSlip'
if (any(prm%B < 0.0_pReal)) extmsg = trim(extmsg)//' B'
if (any(prm%p<=0.0_pReal .or. prm%p>1.0_pReal)) extmsg = trim(extmsg)//' p'
if (any(prm%q< 1.0_pReal .or. prm%q>2.0_pReal)) extmsg = trim(extmsg)//' q'
2018-10-12 20:54:46 +05:30
else slipActive
2019-03-19 03:14:54 +05:30
allocate(prm%b_sl(0))
2018-10-12 20:54:46 +05:30
endif slipActive
2018-06-27 21:08:52 +05:30
2018-10-18 02:43:47 +05:30
!--------------------------------------------------------------------------------------------------
! twin related parameters
2019-03-19 02:54:45 +05:30
prm%N_tw = config%getInts('ntwin', defaultVal=emptyIntArray)
2019-03-20 14:16:49 +05:30
prm%sum_N_tw = sum(prm%N_tw)
if (prm%sum_N_tw > 0) then
2019-03-21 11:59:38 +05:30
prm%P_tw = lattice_SchmidMatrix_twin(prm%N_tw,config%getString('lattice_structure'),&
2019-03-19 02:54:45 +05:30
config%getFloat('c/a',defaultVal=0.0_pReal))
prm%h_tw_tw = lattice_interaction_TwinByTwin(prm%N_tw,&
config%getFloats('interaction_twintwin'), &
config%getString('lattice_structure'))
2018-10-18 02:43:47 +05:30
2019-03-19 03:14:54 +05:30
prm%b_tw = config%getFloats('twinburgers', requiredSize=size(prm%N_tw))
2019-03-19 02:54:45 +05:30
prm%twinsize = config%getFloats('twinsize', requiredSize=size(prm%N_tw))
2019-03-19 03:14:54 +05:30
prm%r = config%getFloats('r_twin', requiredSize=size(prm%N_tw))
2018-10-18 03:26:57 +05:30
2019-03-21 02:21:52 +05:30
prm%xc_twin = config%getFloat('xc_twin')
2019-03-21 11:59:38 +05:30
prm%L_tw = config%getFloat('l0_twin')
2019-03-21 02:21:52 +05:30
prm%i_tw = config%getFloat('cmfptwin')
2018-10-18 03:26:57 +05:30
2019-03-21 11:59:38 +05:30
prm%gamma_char = lattice_characteristicShear_Twin(prm%N_tw,config%getString('lattice_structure'),&
2019-03-19 03:14:54 +05:30
config%getFloat('c/a',defaultVal=0.0_pReal))
2018-12-10 02:50:18 +05:30
2019-03-21 11:59:38 +05:30
prm%C66_tw = lattice_C66_twin(prm%N_tw,prm%C66,config%getString('lattice_structure'),&
2019-03-19 03:14:54 +05:30
config%getFloat('c/a',defaultVal=0.0_pReal))
2018-06-27 21:08:52 +05:30
if (.not. prm%fccTwinTransNucleation) then
prm%Ndot0_twin = config%getFloats('ndot0_twin')
2019-03-19 02:54:45 +05:30
prm%Ndot0_twin = math_expand(prm%Ndot0_twin,prm%N_tw)
endif
2018-10-18 03:26:57 +05:30
! expand: family => system
2019-03-19 03:14:54 +05:30
prm%b_tw = math_expand(prm%b_tw,prm%N_tw)
2019-03-19 02:54:45 +05:30
prm%twinsize = math_expand(prm%twinsize,prm%N_tw)
prm%r = math_expand(prm%r,prm%N_tw)
2018-10-04 18:21:32 +05:30
else
allocate(prm%twinsize(0))
2019-03-19 03:14:54 +05:30
allocate(prm%b_tw(0))
2018-10-04 18:21:32 +05:30
allocate(prm%r(0))
endif
2018-10-04 18:21:32 +05:30
2018-10-19 01:50:26 +05:30
!--------------------------------------------------------------------------------------------------
! transformation related parameters
prm%N_tr = config%getInts('ntrans', defaultVal=emptyIntArray)
2019-03-20 14:16:49 +05:30
prm%sum_N_tr = sum(prm%N_tr)
if (prm%sum_N_tr > 0) then
prm%b_tr = config%getFloats('transburgers')
prm%b_tr = math_expand(prm%b_tr,prm%N_tr)
prm%transStackHeight = config%getFloat('transstackheight', defaultVal=0.0_pReal) ! ToDo: How to handle that???
2019-03-21 02:21:52 +05:30
prm%i_tr = config%getFloat('cmfptrans', defaultVal=0.0_pReal) ! ToDo: How to handle that???
2019-03-21 11:59:38 +05:30
prm%gamma_fcc_hex = config%getFloat('deltag')
prm%xc_trans = config%getFloat('xc_trans', defaultVal=0.0_pReal) ! ToDo: How to handle that???
2019-03-21 11:59:38 +05:30
prm%L_tr = config%getFloat('l0_trans')
2019-03-21 11:59:38 +05:30
prm%h_tr_tr = lattice_interaction_TransByTrans(prm%N_tr,&
config%getFloats('interaction_transtrans'), &
config%getString('lattice_structure'))
2019-03-21 11:59:38 +05:30
prm%C66_tr = lattice_C66_trans(prm%N_tr,prm%C66, &
config%getString('trans_lattice_structure'), &
0.0_pReal, &
config%getFloat('a_bcc', defaultVal=0.0_pReal), &
config%getFloat('a_fcc', defaultVal=0.0_pReal))
2019-03-21 11:59:38 +05:30
prm%P_tr = lattice_SchmidMatrix_trans(prm%N_tr, &
config%getString('trans_lattice_structure'), &
0.0_pReal, &
config%getFloat('a_bcc', defaultVal=0.0_pReal), &
config%getFloat('a_fcc', defaultVal=0.0_pReal))
if (lattice_structure(p) /= LATTICE_fcc_ID) then
prm%Ndot0_trans = config%getFloats('ndot0_trans')
prm%Ndot0_trans = math_expand(prm%Ndot0_trans,prm%N_tr)
endif
2019-01-28 04:06:34 +05:30
prm%lamellarsize = config%getFloats('lamellarsize')
prm%lamellarsize = math_expand(prm%lamellarsize,prm%N_tr)
prm%s = config%getFloats('s_trans',defaultVal=[0.0_pReal])
prm%s = math_expand(prm%s,prm%N_tr)
2018-10-04 18:21:32 +05:30
else
2019-01-28 04:06:34 +05:30
allocate(prm%lamellarsize(0))
2019-03-20 14:16:49 +05:30
allocate(prm%b_tr(0))
endif
2019-03-20 14:16:49 +05:30
if (sum(prm%N_tw) > 0 .or. prm%sum_N_tr > 0) then
prm%SFE_0K = config%getFloat('sfe_0k')
prm%dSFE_dT = config%getFloat('dsfe_dt')
2019-03-21 11:59:38 +05:30
prm%V_cs = config%getFloat('vcrossslip')
endif
2019-03-20 14:16:49 +05:30
if (prm%sum_N_sl > 0 .and. prm%sum_N_tw > 0) then
2019-03-19 02:54:45 +05:30
prm%h_sl_tw = lattice_interaction_SlipByTwin(prm%N_sl,prm%N_tw,&
2019-03-19 02:38:41 +05:30
config%getFloats('interaction_sliptwin'), &
config%getString('lattice_structure'))
2019-03-20 14:16:49 +05:30
if (prm%fccTwinTransNucleation .and. prm%sum_N_tw > 12) write(6,*) 'mist' ! ToDo: implement better test. The model will fail also if N_tw is [6,6]
2018-08-31 19:38:01 +05:30
endif
2019-03-20 14:16:49 +05:30
if (prm%sum_N_sl > 0 .and. prm%sum_N_tr > 0) then
2019-03-21 11:59:38 +05:30
prm%h_sl_tr = lattice_interaction_SlipByTrans(prm%N_sl,prm%N_tr,&
config%getFloats('interaction_sliptrans'), &
config%getString('lattice_structure'))
2019-03-20 14:16:49 +05:30
if (prm%fccTwinTransNucleation .and. prm%sum_N_tr > 12) write(6,*) 'mist' ! ToDo: implement better test. The model will fail also if N_tr is [6,6]
2019-01-29 11:11:27 +05:30
endif
!--------------------------------------------------------------------------------------------------
! shearband related parameters
prm%sbVelocity = config%getFloat('shearbandvelocity',defaultVal=0.0_pReal)
if (prm%sbVelocity > 0.0_pReal) then
prm%sbResistance = config%getFloat('shearbandresistance')
prm%sbQedge = config%getFloat('qedgepersbsystem')
2019-03-21 11:59:38 +05:30
prm%p_sb = config%getFloat('p_shearband')
prm%q_sb = config%getFloat('q_shearband')
2019-01-29 11:11:27 +05:30
! sanity checks
if (prm%sbResistance < 0.0_pReal) extmsg = trim(extmsg)//' shearbandresistance'
if (prm%sbQedge < 0.0_pReal) extmsg = trim(extmsg)//' qedgepersbsystem'
2019-03-21 11:59:38 +05:30
if (prm%p_sb <= 0.0_pReal) extmsg = trim(extmsg)//' p_shearband'
if (prm%q_sb <= 0.0_pReal) extmsg = trim(extmsg)//' q_shearband'
2019-01-29 11:11:27 +05:30
endif
2019-03-21 02:21:52 +05:30
prm%D = config%getFloat('grainsize')
2019-01-29 11:11:27 +05:30
if (config%keyExists('dipoleformationfactor')) call IO_error(1,ext_msg='use /nodipoleformation/')
prm%dipoleformation = .not. config%keyExists('/nodipoleformation/')
2019-01-29 11:11:27 +05:30
2018-10-19 01:04:26 +05:30
!if (Ndot0PerTwinFamily(f,p) < 0.0_pReal) &
2019-03-19 02:47:11 +05:30
! call IO_error(211,el=p,ext_msg='ndot0_twin ('//PLASTICITY_DISLOTWIN_label//')')
2018-10-19 01:04:26 +05:30
if (any(prm%atomicVolume <= 0.0_pReal)) &
2019-03-19 02:47:11 +05:30
call IO_error(211,el=p,ext_msg='cAtomicVolume ('//PLASTICITY_DISLOTWIN_label//')')
2019-03-20 14:16:49 +05:30
if (prm%sum_N_tw > 0) then
2019-03-19 10:55:21 +05:30
if (prm%aTol_rho <= 0.0_pReal) &
call IO_error(211,el=p,ext_msg='aTol_rho ('//PLASTICITY_DISLOTWIN_label//')')
if (prm%aTol_f_tw <= 0.0_pReal) &
call IO_error(211,el=p,ext_msg='aTol_f_tw ('//PLASTICITY_DISLOTWIN_label//')')
2018-10-19 01:04:26 +05:30
endif
2019-03-20 14:16:49 +05:30
if (prm%sum_N_tr > 0) then
2019-03-19 10:55:21 +05:30
if (prm%aTol_f_tr <= 0.0_pReal) &
call IO_error(211,el=p,ext_msg='aTol_f_tr ('//PLASTICITY_DISLOTWIN_label//')')
2018-10-19 01:04:26 +05:30
endif
outputs = config%getStrings('(output)', defaultVal=emptyStringArray)
allocate(prm%outputID(0))
2019-03-19 02:47:11 +05:30
do i= 1, size(outputs)
outputID = undefined_ID
select case(outputs(i))
case ('edge_density')
2019-03-20 14:16:49 +05:30
outputID = merge(rho_mob_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('dipole_density')
2019-03-20 14:16:49 +05:30
outputID = merge(rho_dip_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('shear_rate_slip','shearrate_slip')
2019-03-22 11:18:38 +05:30
outputID = merge(dot_gamma_sl_ID,undefined_ID,prm%sum_N_sl > 0)
2019-03-20 14:16:49 +05:30
outputSize = prm%sum_N_sl
case ('accumulated_shear_slip')
2019-03-20 14:16:49 +05:30
outputID = merge(gamma_sl_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('mfp_slip')
2019-03-20 14:16:49 +05:30
outputID = merge(Lambda_sl_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('resolved_stress_slip')
2019-03-20 14:16:49 +05:30
outputID = merge(resolved_stress_slip_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('threshold_stress_slip')
2019-03-20 14:16:49 +05:30
outputID= merge(threshold_stress_slip_ID,undefined_ID,prm%sum_N_sl > 0)
outputSize = prm%sum_N_sl
case ('twin_fraction')
2019-03-20 14:16:49 +05:30
outputID = merge(f_tw_ID,undefined_ID,prm%sum_N_tw >0)
outputSize = prm%sum_N_tw
case ('mfp_twin')
2019-03-20 14:16:49 +05:30
outputID = merge(Lambda_tw_ID,undefined_ID,prm%sum_N_tw >0)
outputSize = prm%sum_N_tw
case ('resolved_stress_twin')
2019-03-20 14:16:49 +05:30
outputID = merge(resolved_stress_twin_ID,undefined_ID,prm%sum_N_tw >0)
outputSize = prm%sum_N_tw
2019-03-22 13:56:39 +05:30
case ('threshold_stress_twin')
2019-03-22 11:27:32 +05:30
outputID = merge(tau_hat_tw_ID,undefined_ID,prm%sum_N_tw >0)
2019-03-20 14:16:49 +05:30
outputSize = prm%sum_N_tw
case ('strain_trans_fraction')
2019-03-20 14:16:49 +05:30
outputID = f_tr_ID
outputSize = prm%sum_N_tr
end select
if (outputID /= undefined_ID) then
2018-09-05 19:45:57 +05:30
plastic_dislotwin_output(i,phase_plasticityInstance(p)) = outputs(i)
plastic_dislotwin_sizePostResult(i,phase_plasticityInstance(p)) = outputSize
prm%outputID = [prm%outputID, outputID]
endif
2019-01-27 16:07:50 +05:30
enddo
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2019-01-27 16:07:50 +05:30
NipcMyPhase = count(material_phase == p)
2019-03-21 02:21:52 +05:30
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
2019-01-27 16:07:50 +05:30
sizeState = sizeDotState
2019-03-19 02:47:11 +05:30
call material_allocatePlasticState(p,NipcMyPhase,sizeState,sizeDotState,0, &
2019-03-20 14:16:49 +05:30
prm%sum_N_sl,prm%sum_N_tw,prm%sum_N_tr)
2018-09-05 19:45:57 +05:30
plasticState(p)%sizePostResults = sum(plastic_dislotwin_sizePostResult(:,phase_plasticityInstance(p)))
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
! locally defined state aliases and initialization of state0 and aTolState
2019-03-19 02:47:11 +05:30
startIndex = 1
2019-03-20 14:16:49 +05:30
endIndex = prm%sum_N_sl
2019-03-21 02:21:52 +05:30
stt%rho_mob=>plasticState(p)%state(startIndex:endIndex,:)
stt%rho_mob= spread(prm%rho_mob_0,2,NipcMyPhase)
dot%rho_mob=>plasticState(p)%dotState(startIndex:endIndex,:)
2019-03-19 10:55:21 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTol_rho
2018-08-31 15:12:54 +05:30
2019-03-19 02:47:11 +05:30
startIndex = endIndex + 1
2019-03-20 14:16:49 +05:30
endIndex = endIndex + prm%sum_N_sl
2019-03-21 02:21:52 +05:30
stt%rho_dip=>plasticState(p)%state(startIndex:endIndex,:)
stt%rho_dip= spread(prm%rho_dip_0,2,NipcMyPhase)
dot%rho_dip=>plasticState(p)%dotState(startIndex:endIndex,:)
2019-03-19 10:55:21 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTol_rho
2019-01-27 16:07:50 +05:30
2019-03-19 02:47:11 +05:30
startIndex = endIndex + 1
2019-03-20 14:16:49 +05:30
endIndex = endIndex + prm%sum_N_sl
2019-03-21 11:59:38 +05:30
stt%gamma_sl=>plasticState(p)%state(startIndex:endIndex,:)
dot%gamma_sl=>plasticState(p)%dotState(startIndex:endIndex,:)
2019-01-27 16:07:50 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = 1.0e6_pReal !ToDo: better make optional parameter
! global alias
plasticState(p)%slipRate => plasticState(p)%dotState(startIndex:endIndex,:)
plasticState(p)%accumulatedSlip => plasticState(p)%state(startIndex:endIndex,:)
2018-08-31 15:12:54 +05:30
2019-03-19 02:47:11 +05:30
startIndex = endIndex + 1
2019-03-20 14:16:49 +05:30
endIndex = endIndex + prm%sum_N_tw
2019-03-21 02:21:52 +05:30
stt%f_tw=>plasticState(p)%state(startIndex:endIndex,:)
dot%f_tw=>plasticState(p)%dotState(startIndex:endIndex,:)
2019-03-19 10:55:21 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTol_f_tw
2018-08-31 15:12:54 +05:30
2019-03-19 02:47:11 +05:30
startIndex = endIndex + 1
2019-03-20 14:16:49 +05:30
endIndex = endIndex + prm%sum_N_tr
2019-03-21 11:59:38 +05:30
stt%f_tr=>plasticState(p)%state(startIndex:endIndex,:)
dot%f_tr=>plasticState(p)%dotState(startIndex:endIndex,:)
2019-03-19 10:55:21 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTol_f_tr
2019-03-20 14:16:49 +05:30
allocate(dst%Lambda_sl (prm%sum_N_sl, NipcMyPhase),source=0.0_pReal)
allocate(dst%tau_pass (prm%sum_N_sl, NipcMyPhase),source=0.0_pReal)
2018-08-31 15:12:54 +05:30
2019-03-20 14:16:49 +05:30
allocate(dst%Lambda_tw (prm%sum_N_tw, NipcMyPhase),source=0.0_pReal)
2019-03-22 11:27:32 +05:30
allocate(dst%tau_hat_tw (prm%sum_N_tw, NipcMyPhase),source=0.0_pReal)
2019-03-21 02:21:52 +05:30
allocate(dst%tau_r_tw (prm%sum_N_tw, NipcMyPhase),source=0.0_pReal)
allocate(dst%f_tw (prm%sum_N_tw, NipcMyPhase),source=0.0_pReal)
2018-08-31 15:12:54 +05:30
2019-03-20 14:16:49 +05:30
allocate(dst%Lambda_tr (prm%sum_N_tr,NipcMyPhase),source=0.0_pReal)
2019-03-22 11:27:32 +05:30
allocate(dst%tau_hat_tr (prm%sum_N_tr,NipcMyPhase),source=0.0_pReal)
2019-03-21 02:21:52 +05:30
allocate(dst%tau_r_tr (prm%sum_N_tr,NipcMyPhase),source=0.0_pReal)
allocate(dst%f_tr (prm%sum_N_tr,NipcMyPhase),source=0.0_pReal)
2019-01-28 04:06:34 +05:30
2019-01-27 16:07:50 +05:30
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
end associate
2019-01-27 16:07:50 +05:30
enddo
2019-01-27 16:07:50 +05:30
end subroutine plastic_dislotwin_init
2014-07-22 13:13:03 +05:30
2019-01-27 16:07:50 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns the homogenized elasticity matrix
!--------------------------------------------------------------------------------------------------
2019-01-27 16:45:11 +05:30
function plastic_dislotwin_homogenizedC(ipc,ip,el) result(homogenizedC)
use material, only: &
material_phase, &
phase_plasticityInstance, &
phasememberAt
implicit none
real(pReal), dimension(6,6) :: &
2019-01-27 16:45:11 +05:30
homogenizedC
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
2019-03-19 02:47:11 +05:30
integer :: i, &
of
real(pReal) :: f_unrotated
of = phasememberAt(ipc,ip,el)
associate(prm => param(phase_plasticityInstance(material_phase(ipc,ip,el))),&
stt => state(phase_plasticityInstance(material_phase(ipc,ip,el))))
f_unrotated = 1.0_pReal &
2019-03-21 02:21:52 +05:30
- sum(stt%f_tw(1:prm%sum_N_tw,of)) &
2019-03-21 11:59:38 +05:30
- sum(stt%f_tr(1:prm%sum_N_tr,of))
2019-01-27 16:45:11 +05:30
homogenizedC = f_unrotated * prm%C66
2019-03-20 14:16:49 +05:30
do i=1,prm%sum_N_tw
2019-01-27 16:45:11 +05:30
homogenizedC = homogenizedC &
2019-03-21 11:59:38 +05:30
+ stt%f_tw(i,of)*prm%C66_tw(1:6,1:6,i)
enddo
2019-03-20 14:16:49 +05:30
do i=1,prm%sum_N_tr
2019-01-27 16:45:11 +05:30
homogenizedC = homogenizedC &
2019-03-21 11:59:38 +05:30
+ stt%f_tr(i,of)*prm%C66_tr(1:6,1:6,i)
enddo
2019-01-27 16:45:11 +05:30
end associate
2019-01-27 16:45:11 +05:30
end function plastic_dislotwin_homogenizedC
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
subroutine plastic_dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,T,instance,of)
use prec, only: &
2016-05-29 14:15:03 +05:30
tol_math_check, &
dNeq0
use math, only: &
math_eigenValuesVectorsSym, &
2019-03-09 21:28:59 +05:30
math_outer, &
math_symmetric33, &
math_mul33xx33, &
2016-01-09 17:42:31 +05:30
math_mul33x3
implicit none
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
2019-03-19 02:47:11 +05:30
integer, intent(in) :: instance,of
2019-03-21 02:21:52 +05:30
real(pReal), intent(in) :: T
2019-03-19 02:47:11 +05:30
integer :: i,k,l,m,n
real(pReal) :: f_unrotated,StressRatio_p,&
2019-01-27 13:05:07 +05:30
BoltzmannRatio, &
2019-03-20 14:16:49 +05:30
dgamma_dtau, &
2018-08-31 20:06:19 +05:30
tau
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl) :: &
dot_gamma_sl,dgamma_dtau_slip
real(pReal), dimension(param(instance)%sum_N_tw) :: &
dot_gamma_twin,dgamma_dtau_twin
real(pReal), dimension(param(instance)%sum_N_tr) :: &
2019-03-21 02:21:52 +05:30
dot_gamma_tr,dgamma_dtau_trans
2019-03-20 14:16:49 +05:30
real(pReal):: dot_gamma_sb
2019-03-22 11:18:38 +05:30
real(pReal), dimension(3,3) :: eigVectors, P_sb
2019-01-28 00:14:53 +05:30
real(pReal), dimension(3) :: eigValues
logical :: error
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])
2018-09-01 14:15:34 +05:30
associate(prm => param(instance), stt => state(instance))
2018-06-27 21:08:52 +05:30
f_unrotated = 1.0_pReal &
2019-03-21 02:21:52 +05:30
- sum(stt%f_tw(1:prm%sum_N_tw,of)) &
2019-03-21 11:59:38 +05:30
- sum(stt%f_tr(1:prm%sum_N_tr,of))
2014-06-11 17:41:14 +05:30
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
2019-03-21 02:21:52 +05:30
call kinetics_slip(Mp,T,instance,of,dot_gamma_sl,dgamma_dtau_slip)
2019-03-20 14:16:49 +05:30
slipContribution: do i = 1, prm%sum_N_sl
2019-03-21 11:59:38 +05:30
Lp = Lp + dot_gamma_sl(i)*prm%P_sl(1:3,1:3,i)
2019-03-19 02:47:11 +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) &
2019-03-21 11:59:38 +05:30
+ dgamma_dtau_slip(i) * prm%P_sl(k,l,i) * prm%P_sl(m,n,i)
2018-09-01 14:15:34 +05:30
enddo slipContribution
!ToDo: Why do this before shear banding?
Lp = Lp * f_unrotated
dLp_dMp = dLp_dMp * f_unrotated
shearBandingContribution: if(dNeq0(prm%sbVelocity)) then
2018-09-01 14:15:34 +05:30
2019-03-21 02:21:52 +05:30
BoltzmannRatio = prm%sbQedge/(kB*T)
call math_eigenValuesVectorsSym(Mp,eigValues,eigVectors,error)
2018-09-01 14:15:34 +05:30
2019-03-19 02:47:11 +05:30
do i = 1,6
2019-03-22 11:18:38 +05:30
P_sb = 0.5_pReal * math_outer(math_mul33x3(eigVectors,sb_sComposition(1:3,i)),&
2019-03-09 21:28:59 +05:30
math_mul33x3(eigVectors,sb_mComposition(1:3,i)))
2019-03-22 11:18:38 +05:30
tau = math_mul33xx33(Mp,P_sb)
significantShearBandStress: if (abs(tau) > tol_math_check) then
2019-03-21 11:59:38 +05:30
StressRatio_p = (abs(tau)/prm%sbResistance)**prm%p_sb
dot_gamma_sb = sign(prm%sbVelocity*exp(-BoltzmannRatio*(1-StressRatio_p)**prm%q_sb), tau)
dgamma_dtau = abs(dot_gamma_sb)*BoltzmannRatio* prm%p_sb*prm%q_sb/ prm%sbResistance &
* (abs(tau)/prm%sbResistance)**(prm%p_sb-1.0_pReal) &
* (1.0_pReal-StressRatio_p)**(prm%q_sb-1.0_pReal)
2019-03-22 11:18:38 +05:30
Lp = Lp + dot_gamma_sb * P_sb
2019-03-19 02:47:11 +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) &
2019-03-22 11:18:38 +05:30
+ dgamma_dtau * P_sb(k,l) * P_sb(m,n)
endif significantShearBandStress
enddo
2018-09-01 14:15:34 +05:30
endif shearBandingContribution
2019-03-21 02:21:52 +05:30
call kinetics_twin(Mp,T,dot_gamma_sl,instance,of,dot_gamma_twin,dgamma_dtau_twin)
2019-03-20 14:16:49 +05:30
twinContibution: do i = 1, prm%sum_N_tw
2019-03-21 11:59:38 +05:30
Lp = Lp + dot_gamma_twin(i)*prm%P_tw(1:3,1:3,i) * f_unrotated
2019-03-19 02:47:11 +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) &
2019-03-21 11:59:38 +05:30
+ dgamma_dtau_twin(i)* prm%P_tw(k,l,i)*prm%P_tw(m,n,i) * f_unrotated
2018-09-01 14:15:34 +05:30
enddo twinContibution
2019-01-28 00:14:53 +05:30
2019-03-21 02:21:52 +05:30
call kinetics_twin(Mp,T,dot_gamma_sl,instance,of,dot_gamma_tr,dgamma_dtau_trans)
2019-03-20 14:16:49 +05:30
transContibution: do i = 1, prm%sum_N_tr
2019-03-21 11:59:38 +05:30
Lp = Lp + dot_gamma_tr(i)*prm%P_tr(1:3,1:3,i) * f_unrotated
2019-03-19 02:47:11 +05:30
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
2019-01-28 00:14:53 +05:30
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
2019-03-21 11:59:38 +05:30
+ dgamma_dtau_trans(i)* prm%P_tr(k,l,i)*prm%P_tr(m,n,i) * f_unrotated
2019-01-28 00:14:53 +05:30
enddo transContibution
2018-09-01 14:15:34 +05:30
end associate
end subroutine plastic_dislotwin_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief calculates the rate of change of microstructure
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
subroutine plastic_dislotwin_dotState(Mp,T,instance,of)
2013-11-27 21:50:27 +05:30
use prec, only: &
2016-05-29 14:15:03 +05:30
tol_math_check, &
dEq0
2013-11-27 21:50:27 +05:30
use math, only: &
2019-01-29 11:11:27 +05:30
math_clip, &
math_mul33xx33, &
2019-01-29 11:11:27 +05:30
PI
implicit none
real(pReal), dimension(3,3), intent(in):: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T !< temperature at integration point
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
instance, &
of
2019-03-19 02:47:11 +05:30
integer :: i
2019-01-27 13:05:07 +05:30
real(pReal) :: f_unrotated,&
2019-01-29 11:22:55 +05:30
VacancyDiffusion,&
2019-03-22 11:18:38 +05:30
rho_dip_distance, ClimbVelocity, &
2018-08-31 20:06:19 +05:30
tau
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl) :: &
2019-03-22 11:18:38 +05:30
dot_rho_dip_formation, &
dot_rho_dip_climb, &
rho_dip_distance_min, &
2019-03-20 14:16:49 +05:30
dot_gamma_sl
real(pReal), dimension(param(instance)%sum_N_tw) :: &
dot_gamma_twin
real(pReal), dimension(param(instance)%sum_N_tr) :: &
2019-03-21 02:21:52 +05:30
dot_gamma_tr
associate(prm => param(instance), stt => state(instance), &
2019-01-28 02:38:36 +05:30
dot => dotstate(instance), dst => microstructure(instance))
2018-09-01 14:15:34 +05:30
f_unrotated = 1.0_pReal &
2019-03-21 02:21:52 +05:30
- sum(stt%f_tw(1:prm%sum_N_tw,of)) &
2019-03-21 11:59:38 +05:30
- sum(stt%f_tr(1:prm%sum_N_tr,of))
2019-03-21 02:21:52 +05:30
VacancyDiffusion = prm%D0*exp(-prm%Qsd/(kB*T))
2019-03-21 02:21:52 +05:30
call kinetics_slip(Mp,T,instance,of,dot_gamma_sl)
2019-03-21 11:59:38 +05:30
dot%gamma_sl(:,of) = abs(dot_gamma_sl)
2019-01-29 11:22:55 +05:30
2019-03-22 11:18:38 +05:30
rho_dip_distance_min = prm%CEdgeDipMinDistance*prm%b_sl
2019-01-28 00:14:53 +05:30
2019-03-20 14:16:49 +05:30
slipState: do i = 1, prm%sum_N_sl
2019-03-21 11:59:38 +05:30
tau = math_mul33xx33(Mp,prm%P_sl(1:3,1:3,i))
2018-08-31 18:03:42 +05:30
significantSlipStress: if (dEq0(tau)) then
2019-03-22 11:18:38 +05:30
dot_rho_dip_formation(i) = 0.0_pReal
dot_rho_dip_climb(i) = 0.0_pReal
else significantSlipStress
2019-03-22 11:18:38 +05:30
rho_dip_distance = 3.0_pReal*prm%mu*prm%b_sl(i)/(16.0_pReal*PI*abs(tau))
rho_dip_distance = math_clip(rho_dip_distance, right = dst%Lambda_sl(i,of))
rho_dip_distance = math_clip(rho_dip_distance, left = rho_dip_distance_min(i))
if (prm%dipoleFormation) then
2019-03-22 11:18:38 +05:30
dot_rho_dip_formation(i) = 2.0_pReal*(rho_dip_distance-rho_dip_distance_min(i))/prm%b_sl(i) &
* stt%rho_mob(i,of)*abs(dot_gamma_sl(i))
else
2019-03-22 11:18:38 +05:30
dot_rho_dip_formation(i) = 0.0_pReal
endif
2019-03-22 11:18:38 +05:30
if (dEq0(rho_dip_distance-rho_dip_distance_min(i))) then
dot_rho_dip_climb(i) = 0.0_pReal
else
ClimbVelocity = 3.0_pReal*prm%mu*VacancyDiffusion*prm%atomicVolume(i) &
2019-03-22 11:18:38 +05:30
/ (2.0_pReal*PI*kB*T*(rho_dip_distance+rho_dip_distance_min(i)))
dot_rho_dip_climb(i) = 4.0_pReal*ClimbVelocity*stt%rho_dip(i,of) &
/ (rho_dip_distance-rho_dip_distance_min(i))
endif
endif significantSlipStress
enddo slipState
2019-03-22 11:18:38 +05:30
dot%rho_mob(:,of) = abs(dot_gamma_sl)/(prm%b_sl*dst%Lambda_sl(:,of)) &
- dot_rho_dip_formation &
- 2.0_pReal*rho_dip_distance_min/prm%b_sl * stt%rho_mob(:,of)*abs(dot_gamma_sl)
dot%rho_dip(:,of) = dot_rho_dip_formation &
- 2.0_pReal*rho_dip_distance_min/prm%b_sl * stt%rho_dip(:,of)*abs(dot_gamma_sl) &
- dot_rho_dip_climb
2019-03-21 02:21:52 +05:30
call kinetics_twin(Mp,T,dot_gamma_sl,instance,of,dot_gamma_twin)
2019-03-21 11:59:38 +05:30
dot%f_tw(:,of) = f_unrotated*dot_gamma_twin/prm%gamma_char
2018-09-01 14:15:34 +05:30
2019-03-21 02:21:52 +05:30
call kinetics_trans(Mp,T,dot_gamma_sl,instance,of,dot_gamma_tr)
dot%f_tw(:,of) = f_unrotated*dot_gamma_tr
2018-08-31 19:38:01 +05:30
end associate
2019-01-28 00:14:53 +05:30
end subroutine plastic_dislotwin_dotState
2018-09-12 16:55:18 +05:30
2019-01-27 13:05:07 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
subroutine plastic_dislotwin_dependentState(T,instance,of)
2019-01-27 13:05:07 +05:30
use math, only: &
PI
implicit none
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
instance, &
2019-01-27 13:05:07 +05:30
of
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T
2019-01-27 13:05:07 +05:30
2019-03-19 02:47:11 +05:30
integer :: &
2019-01-27 13:05:07 +05:30
i
real(pReal) :: &
sumf_twin,SFE,sumf_trans
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl) :: &
lambda_sl_sl_inv, & !* 1/mean free distance between 2 forest dislocations seen by a moving dislocation
lambda_sl_tw_inv, & !* 1/mean free distance between 2 twin stacks from different systems seen by a moving dislocation
lambda_sl_tr_inv !* 1/mean free distance between 2 martensite lamellar from different systems seen by a moving dislocation
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_tw) :: &
lambda_tw_tw_inv !* 1/mean free distance between 2 twin stacks from different systems seen by a growing twin
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_tr) :: &
lambda_tr_tr_inv !* 1/mean free distance between 2 martensite stacks from different systems seen by a growing martensite (1/lambda_trans)
2019-01-27 13:05:07 +05:30
real(pReal), dimension(:), allocatable :: &
x0, &
fOverStacksize, &
ftransOverLamellarSize
associate(prm => param(instance),&
stt => state(instance),&
2019-01-28 02:38:36 +05:30
dst => microstructure(instance))
2019-01-27 13:05:07 +05:30
2019-03-21 02:21:52 +05:30
sumf_twin = sum(stt%f_tw(1:prm%sum_N_tw,of))
2019-03-21 11:59:38 +05:30
sumf_trans = sum(stt%f_tr(1:prm%sum_N_tr,of))
2019-01-27 13:05:07 +05:30
2019-03-21 02:21:52 +05:30
SFE = prm%SFE_0K + prm%dSFE_dT * T
2019-01-27 13:05:07 +05:30
!* rescaled volume fraction for topology
2019-03-21 02:21:52 +05:30
fOverStacksize = stt%f_tw(1:prm%sum_N_tw,of)/prm%twinsize !ToDo: this is per system
2019-01-28 04:06:34 +05:30
ftransOverLamellarSize = sumf_trans/prm%lamellarsize !ToDo: But this not ...
2019-01-27 13:05:07 +05:30
!Todo: Physically ok, but naming could be adjusted
2019-03-20 14:16:49 +05:30
forall (i = 1:prm%sum_N_sl) &
2019-03-19 19:26:36 +05:30
lambda_sl_sl_inv(i) = &
2019-03-21 02:21:52 +05:30
sqrt(dot_product((stt%rho_mob(1:prm%sum_N_sl,of)+stt%rho_dip(1:prm%sum_N_sl,of)),&
2019-03-20 14:16:49 +05:30
prm%forestProjection(1:prm%sum_N_sl,i)))/prm%CLambdaSlip(i) ! change order and use matmul
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
if (prm%sum_N_tw > 0 .and. prm%sum_N_sl > 0) &
lambda_sl_tw_inv = &
matmul(transpose(prm%h_sl_tw),fOverStacksize)/(1.0_pReal-sumf_twin) ! ToDo: Change order/no transpose
2019-01-27 13:05:07 +05:30
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
!ToDo: needed? if (prm%sum_N_tw > 0) &
lambda_tw_tw_inv = matmul(prm%h_tw_tw,fOverStacksize)/(1.0_pReal-sumf_twin)
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
if (prm%sum_N_tr > 0 .and. prm%sum_N_sl > 0) &
lambda_sl_tr_inv = & ! ToDo: does not work if N_tr is not 12
2019-03-21 11:59:38 +05:30
matmul(transpose(prm%h_sl_tr),ftransOverLamellarSize)/(1.0_pReal-sumf_trans) ! ToDo: remove transpose
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
!ToDo: needed? if (prm%sum_N_tr > 0) &
2019-03-21 11:59:38 +05:30
lambda_tr_tr_inv = matmul(prm%h_tr_tr,ftransOverLamellarSize)/(1.0_pReal-sumf_trans)
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
if ((prm%sum_N_tw > 0) .or. (prm%sum_N_tr > 0)) then ! ToDo: Change order
dst%Lambda_sl(:,of) = &
2019-03-21 02:21:52 +05:30
prm%D/(1.0_pReal+prm%D*&
(lambda_sl_sl_inv + lambda_sl_tw_inv + lambda_sl_tr_inv))
2019-01-27 13:05:07 +05:30
else
2019-03-21 02:21:52 +05:30
dst%Lambda_sl(:,of) = prm%D &
/ (1.0_pReal+prm%D*lambda_sl_sl_inv) !!!!!! correct?
2019-01-27 13:05:07 +05:30
endif
2019-03-21 02:21:52 +05:30
dst%Lambda_tw(:,of) = prm%i_tw*prm%D/(1.0_pReal+prm%D*lambda_tw_tw_inv)
dst%Lambda_tr(:,of) = prm%i_tr*prm%D/(1.0_pReal+prm%D*lambda_tr_tr_inv)
2019-01-27 13:05:07 +05:30
!* threshold stress for dislocation motion
2019-03-20 14:16:49 +05:30
forall (i = 1:prm%sum_N_sl) dst%tau_pass(i,of) = &
2019-03-19 03:14:54 +05:30
prm%mu*prm%b_sl(i)*&
2019-03-21 02:21:52 +05:30
sqrt(dot_product(stt%rho_mob(1:prm%sum_N_sl,of)+stt%rho_dip(1:prm%sum_N_sl,of),&
2019-03-19 02:38:41 +05:30
prm%h_sl_sl(:,i)))
2019-01-27 13:05:07 +05:30
!* threshold stress for growing twin/martensite
2019-03-20 14:16:49 +05:30
if(prm%sum_N_tw == prm%sum_N_sl) &
2019-03-22 11:27:32 +05:30
dst%tau_hat_tw(:,of) = &
2019-03-21 11:59:38 +05:30
(SFE/(3.0_pReal*prm%b_tw)+ 3.0_pReal*prm%b_tw*prm%mu/(prm%L_tw*prm%b_sl)) ! slip burgers here correct?
2019-03-20 14:16:49 +05:30
if(prm%sum_N_tr == prm%sum_N_sl) &
2019-03-22 11:27:32 +05:30
dst%tau_hat_tr(:,of) = &
2019-03-20 14:16:49 +05:30
(SFE/(3.0_pReal*prm%b_tr) + 3.0_pReal*prm%b_tr*prm%mu/&
2019-03-21 11:59:38 +05:30
(prm%L_tr*prm%b_sl) + prm%transStackHeight*prm%gamma_fcc_hex/ (3.0_pReal*prm%b_tr) )
2019-01-27 13:05:07 +05:30
2019-03-21 02:21:52 +05:30
dst%f_tw(:,of) = (PI/4.0_pReal)*prm%twinsize*dst%Lambda_tw(:,of)**2.0_pReal
dst%f_tr(:,of) = (PI/4.0_pReal)*prm%lamellarsize*dst%Lambda_tr(:,of)**2.0_pReal
2019-01-28 04:06:34 +05:30
2019-03-21 11:59:38 +05:30
x0 = prm%mu*prm%b_tw**2.0_pReal/(SFE*8.0_pReal*PI)*(2.0_pReal+prm%nu)/(1.0_pReal-prm%nu) ! ToDo: In the paper, this is the burgers vector for slip and is the same for twin and trans
2019-03-21 02:21:52 +05:30
dst%tau_r_tw(:,of) = prm%mu*prm%b_tw/(2.0_pReal*PI)*(1.0_pReal/(x0+prm%xc_twin)+cos(pi/3.0_pReal)/x0)
2019-01-27 13:05:07 +05:30
2019-03-20 14:16:49 +05:30
x0 = prm%mu*prm%b_tr**2.0_pReal/(SFE*8.0_pReal*PI)*(2.0_pReal+prm%nu)/(1.0_pReal-prm%nu) ! ToDo: In the paper, this is the burgers vector for slip
2019-03-21 02:21:52 +05:30
dst%tau_r_tr(:,of) = prm%mu*prm%b_tr/(2.0_pReal*PI)*(1.0_pReal/(x0+prm%xc_trans)+cos(pi/3.0_pReal)/x0)
2019-01-27 13:05:07 +05:30
2019-01-28 00:14:53 +05:30
end associate
2019-01-27 13:05:07 +05:30
end subroutine plastic_dislotwin_dependentState
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
function plastic_dislotwin_postResults(Mp,T,instance,of) result(postResults)
use prec, only: &
tol_math_check, &
dEq0
use math, only: &
PI, &
2019-01-27 16:07:50 +05:30
math_mul33xx33
implicit none
real(pReal), dimension(3,3),intent(in) :: &
Mp !< 2nd Piola Kirchhoff stress tensor in Mandel notation
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T !< temperature at integration point
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
instance, &
of
real(pReal), dimension(sum(plastic_dislotwin_sizePostResult(:,instance))) :: &
postResults
2019-03-19 02:47:11 +05:30
integer :: &
2019-01-28 00:14:53 +05:30
o,c,j
2019-01-28 02:38:36 +05:30
associate(prm => param(instance), stt => state(instance), dst => microstructure(instance))
2019-03-19 02:47:11 +05:30
c = 0
2019-01-28 00:14:53 +05:30
2019-03-19 02:47:11 +05:30
do o = 1,size(prm%outputID)
select case(prm%outputID(o))
2019-03-19 02:38:41 +05:30
case (rho_mob_ID)
2019-03-21 02:21:52 +05:30
postResults(c+1:c+prm%sum_N_sl) = stt%rho_mob(1:prm%sum_N_sl,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_sl
2019-03-19 02:38:41 +05:30
case (rho_dip_ID)
2019-03-21 02:21:52 +05:30
postResults(c+1:c+prm%sum_N_sl) = stt%rho_dip(1:prm%sum_N_sl,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_sl
2019-03-22 11:18:38 +05:30
case (dot_gamma_sl_ID)
2019-03-21 02:21:52 +05:30
call kinetics_slip(Mp,T,instance,of,postResults(c+1:c+prm%sum_N_sl))
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_sl
2019-03-19 02:38:41 +05:30
case (gamma_sl_ID)
2019-03-21 11:59:38 +05:30
postResults(c+1:c+prm%sum_N_sl) = stt%gamma_sl(1:prm%sum_N_sl,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_sl
case (Lambda_sl_ID)
postResults(c+1:c+prm%sum_N_sl) = dst%Lambda_sl(1:prm%sum_N_sl,of)
c = c + prm%sum_N_sl
case (resolved_stress_slip_ID)
2019-03-20 14:16:49 +05:30
do j = 1, prm%sum_N_sl
2019-03-21 11:59:38 +05:30
postResults(c+j) = math_mul33xx33(Mp,prm%P_sl(1:3,1:3,j))
enddo
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_sl
case (threshold_stress_slip_ID)
2019-03-20 14:16:49 +05:30
postResults(c+1:c+prm%sum_N_sl) = dst%tau_pass(1:prm%sum_N_sl,of)
c = c + prm%sum_N_sl
2019-03-19 02:38:41 +05:30
case (f_tw_ID)
2019-03-21 02:21:52 +05:30
postResults(c+1:c+prm%sum_N_tw) = stt%f_tw(1:prm%sum_N_tw,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_tw
case (Lambda_tw_ID)
postResults(c+1:c+prm%sum_N_tw) = dst%Lambda_tw(1:prm%sum_N_tw,of)
c = c + prm%sum_N_tw
case (resolved_stress_twin_ID)
2019-03-20 14:16:49 +05:30
do j = 1, prm%sum_N_tw
2019-03-21 11:59:38 +05:30
postResults(c+j) = math_mul33xx33(Mp,prm%P_tw(1:3,1:3,j))
enddo
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_tw
2019-03-22 11:27:32 +05:30
case (tau_hat_tw_ID)
postResults(c+1:c+prm%sum_N_tw) = dst%tau_hat_tw(1:prm%sum_N_tw,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_tw
2019-03-20 14:16:49 +05:30
case (f_tr_ID)
2019-03-21 11:59:38 +05:30
postResults(c+1:c+prm%sum_N_tr) = stt%f_tr(1:prm%sum_N_tr,of)
2019-03-20 14:16:49 +05:30
c = c + prm%sum_N_tr
end select
enddo
2019-01-28 00:14:53 +05:30
end associate
2019-01-28 00:14:53 +05:30
end function plastic_dislotwin_postResults
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine plastic_dislotwin_results(instance,group)
#if defined(PETSc) || defined(DAMASKHDF5)
use results
implicit none
integer, intent(in) :: instance
character(len=*) :: group
integer :: o
associate(prm => param(instance), stt => state(instance))
2019-03-19 02:47:11 +05:30
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
end select
enddo outputsLoop
end associate
#else
integer, intent(in) :: instance
character(len=*) :: group
#endif
end subroutine plastic_dislotwin_results
2018-09-12 16:55:18 +05:30
!--------------------------------------------------------------------------------------------------
2019-01-27 21:29:44 +05:30
!> @brief Shear rates on slip systems, their derivatives with respect to resolved stress and the
! resolved stresss
!> @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
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
pure subroutine kinetics_slip(Mp,T,instance,of, &
2019-03-20 14:16:49 +05:30
dot_gamma_sl,dgamma_dtau_slip,tau_slip)
2018-09-12 16:55:18 +05:30
use prec, only: &
tol_math_check, &
dNeq0
use math, only: &
math_mul33xx33
implicit none
2019-01-27 21:29:44 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T !< temperature
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
2019-01-27 21:29:44 +05:30
instance, &
2018-09-12 16:55:18 +05:30
of
2019-01-27 21:29:44 +05:30
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl), intent(out) :: &
dot_gamma_sl
real(pReal), dimension(param(instance)%sum_N_sl), optional, intent(out) :: &
dgamma_dtau_slip, &
2019-01-27 21:29:44 +05:30
tau_slip
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl) :: &
dgamma_dtau
2018-09-12 16:55:18 +05:30
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl) :: &
tau, &
2018-09-12 16:55:18 +05:30
stressRatio, &
StressRatio_p, &
BoltzmannRatio, &
2019-01-27 21:29:44 +05:30
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, &
2019-01-27 21:29:44 +05:30
tau_eff !< effective resolved stress
2019-03-19 02:47:11 +05:30
integer :: i
2019-01-27 21:29:44 +05:30
associate(prm => param(instance), stt => state(instance), dst => microstructure(instance))
2018-09-12 16:55:18 +05:30
2019-03-20 14:16:49 +05:30
do i = 1, prm%sum_N_sl
2019-03-21 11:59:38 +05:30
tau(i) = math_mul33xx33(Mp,prm%P_sl(1:3,1:3,i))
2018-09-12 16:55:18 +05:30
enddo
tau_eff = abs(tau)-dst%tau_pass(:,of)
significantStress: where(tau_eff > tol_math_check)
2019-03-22 11:27:32 +05:30
stressRatio = tau_eff/prm%tau_0
StressRatio_p = stressRatio** prm%p
2019-03-21 02:21:52 +05:30
BoltzmannRatio = prm%Delta_F/(kB*T)
v_wait_inverse = prm%v0**(-1.0_pReal) * exp(BoltzmannRatio*(1.0_pReal-StressRatio_p)** prm%q)
2019-03-19 03:14:54 +05:30
v_run_inverse = prm%B/(tau_eff*prm%b_sl)
2019-03-21 02:21:52 +05:30
dot_gamma_sl = sign(stt%rho_mob(:,of)*prm%b_sl/(v_wait_inverse+v_run_inverse),tau)
2018-09-12 16:55:18 +05:30
dV_wait_inverse_dTau = v_wait_inverse * prm%p * prm%q * BoltzmannRatio &
* (stressRatio**(prm%p-1.0_pReal)) &
* (1.0_pReal-StressRatio_p)**(prm%q-1.0_pReal) &
2019-03-22 11:27:32 +05:30
/ prm%tau_0
dV_run_inverse_dTau = v_run_inverse/tau_eff
dV_dTau = (dV_wait_inverse_dTau+dV_run_inverse_dTau) &
/ (v_wait_inverse+v_run_inverse)**2.0_pReal
2019-03-21 02:21:52 +05:30
dgamma_dtau = dV_dTau*stt%rho_mob(:,of)*prm%b_sl
else where significantStress
2019-03-20 14:16:49 +05:30
dot_gamma_sl = 0.0_pReal
dgamma_dtau = 0.0_pReal
end where significantStress
2019-01-28 00:14:53 +05:30
end associate
2019-03-20 14:16:49 +05:30
if(present(dgamma_dtau_slip)) dgamma_dtau_slip = dgamma_dtau
2019-01-27 21:29:44 +05:30
if(present(tau_slip)) tau_slip = tau
end subroutine kinetics_slip
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on twin systems
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
pure subroutine kinetics_twin(Mp,T,dot_gamma_sl,instance,of,&
2019-03-20 14:16:49 +05:30
dot_gamma_twin,dgamma_dtau_twin)
use prec, only: &
tol_math_check, &
dNeq0
use math, only: &
math_mul33xx33
implicit none
2019-01-28 00:14:53 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T !< temperature
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
2019-01-28 00:14:53 +05:30
instance, &
of
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl), intent(in) :: &
dot_gamma_sl
2019-01-28 00:14:53 +05:30
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_tw), intent(out) :: &
dot_gamma_twin
real(pReal), dimension(param(instance)%sum_N_tw), optional, intent(out) :: &
dgamma_dtau_twin
2018-09-12 16:55:18 +05:30
2019-03-20 14:16:49 +05:30
real, dimension(param(instance)%sum_N_tw) :: &
tau, &
2019-01-28 00:14:53 +05:30
Ndot0, &
stressRatio_r, &
2019-03-20 14:16:49 +05:30
dgamma_dtau
2019-03-19 02:47:11 +05:30
integer :: i,s1,s2
2019-01-28 00:14:53 +05:30
associate(prm => param(instance), stt => state(instance), dst => microstructure(instance))
2019-03-20 14:16:49 +05:30
do i = 1, prm%sum_N_tw
2019-03-21 11:59:38 +05:30
tau(i) = math_mul33xx33(Mp,prm%P_tw(1:3,1:3,i))
isFCC: if (prm%fccTwinTransNucleation) then
s1=prm%fcc_twinNucleationSlipPair(1,i)
s2=prm%fcc_twinNucleationSlipPair(2,i)
2019-03-21 02:21:52 +05:30
if (tau(i) < dst%tau_r_tw(i,of)) then
Ndot0=(abs(dot_gamma_sl(s1))*(stt%rho_mob(s2,of)+stt%rho_dip(s2,of))+&
abs(dot_gamma_sl(s2))*(stt%rho_mob(s1,of)+stt%rho_dip(s1,of)))/& ! ToDo: MD: it would be more consistent to use shearrates from state
2019-03-21 11:59:38 +05:30
(prm%L_tw*prm%b_sl(i))*&
(1.0_pReal-exp(-prm%V_cs/(kB*T)*&
2019-03-21 02:21:52 +05:30
(dst%tau_r_tw(i,of)-tau)))
else
2019-01-28 00:14:53 +05:30
Ndot0=0.0_pReal
end if
else isFCC
2019-01-28 00:14:53 +05:30
Ndot0=prm%Ndot0_twin(i)
endif isFCC
enddo
significantStress: where(tau > tol_math_check)
2019-03-22 11:27:32 +05:30
StressRatio_r = (dst%tau_hat_tw(:,of)/tau)**prm%r
2019-03-21 11:59:38 +05:30
dot_gamma_twin = prm%gamma_char * dst%f_tw(:,of) * Ndot0*exp(-StressRatio_r)
2019-03-20 14:16:49 +05:30
dgamma_dtau = (dot_gamma_twin*prm%r/tau)*StressRatio_r
else where significantStress
2019-03-20 14:16:49 +05:30
dot_gamma_twin = 0.0_pReal
dgamma_dtau = 0.0_pReal
end where significantStress
2019-01-28 00:14:53 +05:30
end associate
2019-03-20 14:16:49 +05:30
if(present(dgamma_dtau_twin)) dgamma_dtau_twin = dgamma_dtau
2018-09-12 16:55:18 +05:30
end subroutine kinetics_twin
2019-01-28 00:14:53 +05:30
2018-09-15 14:13:05 +05:30
!--------------------------------------------------------------------------------------------------
2019-01-28 00:14:53 +05:30
!> @brief calculates shear rates on twin systems
2018-09-15 14:13:05 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-21 02:21:52 +05:30
pure subroutine kinetics_trans(Mp,T,dot_gamma_sl,instance,of,&
dot_gamma_tr,dgamma_dtau_trans)
2018-09-15 14:13:05 +05:30
use prec, only: &
tol_math_check, &
dNeq0
use math, only: &
math_mul33xx33
implicit none
2019-01-28 00:14:53 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
2019-03-21 02:21:52 +05:30
T !< temperature
2019-03-19 02:47:11 +05:30
integer, intent(in) :: &
2019-01-28 00:14:53 +05:30
instance, &
2018-09-15 14:13:05 +05:30
of
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_sl), intent(in) :: &
dot_gamma_sl
2019-01-28 00:14:53 +05:30
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_tr), intent(out) :: &
2019-03-21 02:21:52 +05:30
dot_gamma_tr
2019-03-20 14:16:49 +05:30
real(pReal), dimension(param(instance)%sum_N_tr), optional, intent(out) :: &
dgamma_dtau_trans
2018-09-15 14:13:05 +05:30
2019-03-20 14:16:49 +05:30
real, dimension(param(instance)%sum_N_tr) :: &
2018-09-15 14:13:05 +05:30
tau, &
2019-01-28 00:14:53 +05:30
Ndot0, &
stressRatio_s, &
2019-03-20 14:16:49 +05:30
dgamma_dtau
2018-09-15 14:13:05 +05:30
2019-03-19 02:47:11 +05:30
integer :: i,s1,s2
2019-01-28 00:14:53 +05:30
associate(prm => param(instance), stt => state(instance), dst => microstructure(instance))
2018-09-15 14:13:05 +05:30
2019-03-20 14:16:49 +05:30
do i = 1, prm%sum_N_tr
2019-03-21 11:59:38 +05:30
tau(i) = math_mul33xx33(Mp,prm%P_tr(1:3,1:3,i))
isFCC: if (prm%fccTwinTransNucleation) then
2018-09-15 14:13:05 +05:30
s1=prm%fcc_twinNucleationSlipPair(1,i)
s2=prm%fcc_twinNucleationSlipPair(2,i)
2019-03-21 02:21:52 +05:30
if (tau(i) < dst%tau_r_tr(i,of)) then
Ndot0=(abs(dot_gamma_sl(s1))*(stt%rho_mob(s2,of)+stt%rho_dip(s2,of))+&
abs(dot_gamma_sl(s2))*(stt%rho_mob(s1,of)+stt%rho_dip(s1,of)))/& ! ToDo: MD: it would be more consistent to use shearrates from state
2019-03-21 11:59:38 +05:30
(prm%L_tr*prm%b_sl(i))*&
(1.0_pReal-exp(-prm%V_cs/(kB*T)*&
2019-03-21 02:21:52 +05:30
(dst%tau_r_tr(i,of)-tau)))
2018-09-15 14:13:05 +05:30
else
2019-01-28 00:14:53 +05:30
Ndot0=0.0_pReal
2018-09-15 14:13:05 +05:30
end if
else isFCC
2019-01-28 00:14:53 +05:30
Ndot0=prm%Ndot0_trans(i)
2018-09-15 14:13:05 +05:30
endif isFCC
enddo
2019-01-28 00:14:53 +05:30
significantStress: where(tau > tol_math_check)
2019-03-22 11:27:32 +05:30
StressRatio_s = (dst%tau_hat_tr(:,of)/tau)**prm%s
2019-03-21 02:21:52 +05:30
dot_gamma_tr = dst%f_tr(:,of) * Ndot0*exp(-StressRatio_s)
dgamma_dtau = (dot_gamma_tr*prm%r/tau)*StressRatio_s
2019-01-28 00:14:53 +05:30
else where significantStress
2019-03-21 02:21:52 +05:30
dot_gamma_tr = 0.0_pReal
2019-03-20 14:16:49 +05:30
dgamma_dtau = 0.0_pReal
2019-01-28 00:14:53 +05:30
end where significantStress
end associate
2019-03-20 14:16:49 +05:30
if(present(dgamma_dtau_trans)) dgamma_dtau_trans = dgamma_dtau
2019-01-28 00:14:53 +05:30
end subroutine kinetics_trans
end module plastic_dislotwin