DAMASK_EICMD/src/plastic_disloUCLA.f90

709 lines
31 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @author David Cereceda, Lawrence Livermore National Laboratory
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief crystal plasticity model for bcc metals, especially Tungsten
!--------------------------------------------------------------------------------------------------
2015-01-21 20:44:00 +05:30
module plastic_disloUCLA
use prec, only: &
2019-03-19 03:27:23 +05:30
pReal
implicit none
private
2019-03-19 03:23:24 +05:30
integer, dimension(:,:), allocatable, target, public :: &
2016-01-16 12:36:34 +05:30
plastic_disloUCLA_sizePostResult !< size of each post result output
2019-03-19 03:27:23 +05:30
character(len=64), dimension(:,:), allocatable, target, public :: &
2016-01-16 12:36:34 +05:30
plastic_disloUCLA_output !< name of each post result output
real(pReal), parameter, private :: &
2016-04-25 02:06:35 +05:30
kB = 1.38e-23_pReal !< Boltzmann constant in J/Kelvin
2018-12-05 03:00:07 +05:30
enum, bind(c)
enumerator :: &
undefined_ID, &
2019-03-19 10:16:47 +05:30
rho_mob_ID, &
rho_dip_ID, &
gamma_dot_sl_ID, &
gamma_sl_ID, &
2019-03-19 10:40:23 +05:30
Lambda_sl_ID, &
2019-01-28 02:43:45 +05:30
thresholdstress_ID
end enum
2018-11-28 00:30:45 +05:30
type, private :: tParameters
real(pReal) :: &
2018-11-29 13:14:31 +05:30
aTolRho, &
D, & !< grain size
2018-12-03 15:55:29 +05:30
mu, &
2018-12-21 10:45:01 +05:30
D0, & !< prefactor for self-diffusion coefficient
Qsd !< activation energy for dislocation climb
2019-03-19 03:23:24 +05:30
real(pReal), dimension(:), allocatable :: &
rho_mob_0, & !< initial dislocation density
rho_dip_0, & !< initial dipole density
b_sl, & !< magnitude of burgers vector [m]
2018-12-09 18:59:19 +05:30
nonSchmidCoeff, &
minDipDistance, &
2018-12-21 20:31:16 +05:30
CLambda, & !< Adj. parameter for distance between 2 forest dislocations
2018-12-09 18:59:19 +05:30
atomicVolume, &
2019-03-19 10:40:23 +05:30
tau_0, &
2018-12-09 18:59:19 +05:30
!* mobility law parameters
2019-03-19 10:40:23 +05:30
delta_F, & !< activation energy for glide [J]
2018-12-21 20:31:16 +05:30
v0, & !< dislocation velocity prefactor [m/s]
2018-12-21 10:45:01 +05:30
p, & !< p-exponent in glide velocity
q, & !< q-exponent in glide velocity
B, & !< friction coefficient
kink_height, & !< height of the kink pair
2018-12-21 20:31:16 +05:30
w, & !< width of the kink pair
omega !< attempt frequency for kink pair nucleation
2019-03-19 03:23:24 +05:30
real(pReal), dimension(:,:), allocatable :: &
2019-03-19 10:40:23 +05:30
h_sl_sl, & !< slip resistance from slip activity
2018-12-09 18:59:19 +05:30
forestProjectionEdge
2019-03-19 03:27:23 +05:30
real(pReal), dimension(:,:,:), allocatable :: &
Schmid, &
2018-11-28 00:30:45 +05:30
nonSchmid_pos, &
nonSchmid_neg
2019-03-19 03:23:24 +05:30
integer :: &
sum_N_sl !< total number of active slip system
integer, dimension(:), allocatable :: &
2019-03-19 10:40:23 +05:30
N_sl !< number of active slip systems for each family
2019-03-19 03:23:24 +05:30
integer(kind(undefined_ID)), dimension(:),allocatable :: &
2018-11-28 00:30:45 +05:30
outputID !< ID of each post result output
2018-12-04 04:36:46 +05:30
logical :: &
2019-01-27 16:07:50 +05:30
dipoleFormation !< flag indicating consideration of dipole formation
2018-11-28 00:30:45 +05:30
end type !< container type for internal constitutive parameters
2018-12-05 03:00:07 +05:30
type, private :: tDisloUCLAState
2019-03-19 03:23:24 +05:30
real(pReal), dimension(:,:), pointer :: &
rho_mob, &
rho_dip, &
2019-03-19 03:27:23 +05:30
gamma_sl
end type tDisloUCLAState
2018-12-05 03:00:07 +05:30
type, private :: tDisloUCLAdependentState
2019-03-19 03:23:24 +05:30
real(pReal), dimension(:,:), allocatable :: &
mfp, &
2018-12-21 17:16:43 +05:30
dislocationSpacing, &
threshold_stress
2018-12-05 03:00:07 +05:30
end type tDisloUCLAdependentState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:), private :: param
type(tDisloUCLAState), allocatable, dimension(:), private :: &
dotState, &
state
type(tDisloUCLAdependentState), allocatable, dimension(:), private :: dependentState
2018-12-05 03:00:07 +05:30
public :: &
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_init, &
plastic_disloUCLA_dependentState, &
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_LpAndItsTangent, &
plastic_disloUCLA_dotState, &
plastic_disloUCLA_postResults, &
plastic_disloUCLA_results
private :: &
2018-11-29 13:51:58 +05:30
kinetics
contains
2019-01-07 12:34:02 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2018-12-05 02:03:32 +05:30
subroutine plastic_disloUCLA_init()
2018-12-09 22:23:20 +05:30
use prec, only: &
pStringLen
use debug, only: &
debug_level,&
debug_constitutive,&
debug_levelBasic
use math, only: &
math_expand
use IO, only: &
2019-03-09 05:37:26 +05:30
IO_error
use material, only: &
phase_plasticity, &
phase_plasticityInstance, &
phase_Noutput, &
material_allocatePlasticState, &
2015-01-15 16:26:15 +05:30
PLASTICITY_DISLOUCLA_label, &
PLASTICITY_DISLOUCLA_ID, &
material_phase, &
plasticState
use config, only: &
2018-11-28 00:30:45 +05:30
config_phase
use lattice
2018-12-05 03:00:07 +05:30
implicit none
2019-03-19 03:23:24 +05:30
integer :: &
2018-12-21 10:45:01 +05:30
Ninstance, &
p, i, &
2019-01-07 12:34:02 +05:30
NipcMyPhase, &
2018-12-21 10:45:01 +05:30
sizeState, sizeDotState, &
startIndex, endIndex
2019-03-19 03:27:23 +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-12-05 03:00:07 +05:30
integer(kind(undefined_ID)) :: &
outputID
character(len=pStringLen) :: &
extmsg = ''
character(len=65536), dimension(:), allocatable :: &
outputs
write(6,'(/,a)') ' <<<+- plastic_'//PLASTICITY_DISLOUCLA_label//' init -+>>>'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Cereceda et al., International Journal of Plasticity 78:242256, 2016'
write(6,'(a)') ' https://dx.doi.org/10.1016/j.ijplas.2015.09.002'
2018-12-05 03:00:07 +05:30
2019-03-19 03:27:23 +05:30
Ninstance = count(phase_plasticity == PLASTICITY_DISLOUCLA_ID)
2019-03-19 03:23:24 +05:30
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
2018-12-09 22:06:01 +05:30
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-03-19 03:23:24 +05:30
allocate(plastic_disloUCLA_sizePostResult(maxval(phase_Noutput),Ninstance),source=0)
2018-12-09 22:06:01 +05:30
allocate(plastic_disloUCLA_output(maxval(phase_Noutput),Ninstance))
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_output = ''
2018-12-03 15:55:29 +05:30
2018-12-09 22:06:01 +05:30
allocate(param(Ninstance))
allocate(state(Ninstance))
allocate(dotState(Ninstance))
allocate(dependentState(Ninstance))
2018-11-28 00:30:45 +05:30
2019-03-19 03:23:24 +05:30
do p = 1, size(phase_plasticity)
2018-11-28 00:30:45 +05:30
if (phase_plasticity(p) /= PLASTICITY_DISLOUCLA_ID) cycle
associate(prm => param(phase_plasticityInstance(p)), &
dot => dotState(phase_plasticityInstance(p)), &
2018-12-04 04:36:46 +05:30
stt => state(phase_plasticityInstance(p)), &
dst => dependentState(phase_plasticityInstance(p)), &
config => config_phase(p))
2018-11-28 00:30:45 +05:30
!--------------------------------------------------------------------------------------------------
! optional parameters that need to be defined
2018-11-30 12:55:23 +05:30
prm%mu = lattice_mu(p)
2018-11-28 00:30:45 +05:30
prm%aTolRho = config%getFloat('atol_rho')
! sanity checks
if (prm%aTolRho <= 0.0_pReal) extmsg = trim(extmsg)//' atol_rho'
2018-11-28 00:30:45 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
2019-03-19 10:40:23 +05:30
prm%N_sl = config%getInts('nslip',defaultVal=emptyIntArray)
prm%sum_N_sl = sum(prm%N_sl)
slipActive: if (prm%sum_N_sl > 0) then
2019-03-19 10:40:23 +05:30
prm%Schmid = lattice_SchmidMatrix_slip(prm%N_sl,config%getString('lattice_structure'),&
config%getFloat('c/a',defaultVal=0.0_pReal))
if(trim(config%getString('lattice_structure')) == 'bcc') then
prm%nonSchmidCoeff = config%getFloats('nonschmid_coefficients',&
2019-01-07 12:34:02 +05:30
defaultVal = emptyRealArray)
2019-03-19 10:40:23 +05:30
prm%nonSchmid_pos = lattice_nonSchmidMatrix(prm%N_sl,prm%nonSchmidCoeff,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(prm%N_sl,prm%nonSchmidCoeff,-1)
2018-11-28 00:30:45 +05:30
else
prm%nonSchmid_pos = prm%Schmid
prm%nonSchmid_neg = prm%Schmid
2018-11-28 00:30:45 +05:30
endif
2019-01-27 16:07:50 +05:30
2019-03-19 10:40:23 +05:30
prm%h_sl_sl = lattice_interaction_SlipBySlip(prm%N_sl, &
config%getFloats('interaction_slipslip'), &
config%getString('lattice_structure'))
prm%forestProjectionEdge = 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))
2019-03-19 10:40:23 +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 10:40:23 +05:30
prm%v0 = config%getFloats('v0', requiredSize=size(prm%N_sl))
prm%b_sl = config%getFloats('slipburgers', requiredSize=size(prm%N_sl))
2019-03-19 10:40:23 +05:30
prm%delta_F = config%getFloats('qedge', requiredSize=size(prm%N_sl))
prm%clambda = config%getFloats('clambdaslip', requiredSize=size(prm%N_sl))
prm%tau_0 = config%getFloats('tau_peierls', requiredSize=size(prm%N_sl))
prm%p = config%getFloats('p_slip', requiredSize=size(prm%N_sl), &
defaultVal=[(1.0_pReal,i=1,size(prm%N_sl))])
prm%q = config%getFloats('q_slip', requiredSize=size(prm%N_sl), &
defaultVal=[(1.0_pReal,i=1,size(prm%N_sl))])
prm%kink_height = config%getFloats('kink_height', requiredSize=size(prm%N_sl))
prm%w = config%getFloats('kink_width', requiredSize=size(prm%N_sl))
prm%omega = config%getFloats('omega', requiredSize=size(prm%N_sl))
prm%B = config%getFloats('friction_coeff', requiredSize=size(prm%N_sl))
prm%D = config%getFloat('grainsize')
prm%D0 = config%getFloat('d0')
prm%Qsd = config%getFloat('qsd')
prm%atomicVolume = config%getFloat('catomicvolume') * prm%b_sl**3.0_pReal
prm%minDipDistance = config%getFloat('cedgedipmindistance') * prm%b_sl
prm%dipoleformation = config%getFloat('dipoleformationfactor') > 0.0_pReal !should be on by default, ToDo: change to /key/-type key
2018-11-29 15:01:02 +05:30
! expand: family => system
2019-03-19 10:40:23 +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 10:40:23 +05:30
prm%q = math_expand(prm%q, prm%N_sl)
prm%p = math_expand(prm%p, prm%N_sl)
prm%delta_F = math_expand(prm%delta_F, prm%N_sl)
prm%b_sl = math_expand(prm%b_sl, prm%N_sl)
2019-03-19 10:40:23 +05:30
prm%kink_height = math_expand(prm%kink_height, prm%N_sl)
prm%w = math_expand(prm%w, prm%N_sl)
prm%omega = math_expand(prm%omega, prm%N_sl)
prm%tau_0 = math_expand(prm%tau_0, prm%N_sl)
prm%v0 = math_expand(prm%v0, prm%N_sl)
prm%B = math_expand(prm%B, prm%N_sl)
prm%clambda = math_expand(prm%clambda, prm%N_sl)
prm%atomicVolume = math_expand(prm%atomicVolume, prm%N_sl)
prm%minDipDistance = math_expand(prm%minDipDistance, prm%N_sl)
2018-12-09 22:05:48 +05:30
2018-12-09 22:23:20 +05:30
! sanity checks
2019-01-07 12:34:02 +05:30
if ( prm%D0 <= 0.0_pReal) extmsg = trim(extmsg)//' d0'
if ( prm%Qsd <= 0.0_pReal) extmsg = trim(extmsg)//' qsd'
if (any(prm%rho_mob_0 < 0.0_pReal)) extmsg = trim(extmsg)//' rhoedge0'
if (any(prm%rho_dip_0 < 0.0_pReal)) extmsg = trim(extmsg)//' rhoedgedip0'
2019-01-07 12:34:02 +05:30
if (any(prm%v0 < 0.0_pReal)) extmsg = trim(extmsg)//' v0'
if (any(prm%b_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' slipb_sl'
2019-03-19 10:40:23 +05:30
if (any(prm%delta_F <= 0.0_pReal)) extmsg = trim(extmsg)//' qedge'
if (any(prm%tau_0 < 0.0_pReal)) extmsg = trim(extmsg)//' tau_0'
if (any(prm%minDipDistance <= 0.0_pReal)) extmsg = trim(extmsg)//' cedgedipmindistance or slipb_sl'
if (any(prm%atomicVolume <= 0.0_pReal)) extmsg = trim(extmsg)//' catomicvolume or slipb_sl'
2019-01-07 12:34:02 +05:30
else slipActive
allocate(prm%rho_mob_0(0))
allocate(prm%rho_dip_0(0))
2018-11-28 00:30:45 +05:30
endif slipActive
2018-11-28 10:29:03 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') &
2019-03-19 03:23:24 +05:30
call IO_error(211,ext_msg=trim(extmsg)//'('//PLASTICITY_DISLOUCLA_label//')')
2018-11-28 10:29:03 +05:30
!--------------------------------------------------------------------------------------------------
! output pararameters
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
2018-11-28 10:29:03 +05:30
allocate(prm%outputID(0))
2019-03-19 03:23:24 +05:30
do i=1, size(outputs)
2018-11-28 10:29:03 +05:30
outputID = undefined_ID
select case(trim(outputs(i)))
2018-12-05 03:00:07 +05:30
case ('edge_density')
outputID = merge(rho_mob_ID,undefined_ID,prm%sum_N_sl>0)
2018-12-05 03:00:07 +05:30
case ('dipole_density')
outputID = merge(rho_dip_ID,undefined_ID,prm%sum_N_sl>0)
2018-12-05 03:00:07 +05:30
case ('shear_rate','shearrate','shear_rate_slip','shearrate_slip')
outputID = merge(gamma_dot_sl_ID,undefined_ID,prm%sum_N_sl>0)
2018-12-05 03:00:07 +05:30
case ('accumulated_shear','accumulatedshear','accumulated_shear_slip')
outputID = merge(gamma_sl_ID,undefined_ID,prm%sum_N_sl>0)
2018-12-05 03:00:07 +05:30
case ('mfp','mfp_slip')
outputID = merge(Lambda_sl_ID,undefined_ID,prm%sum_N_sl>0)
2018-12-05 03:00:07 +05:30
case ('threshold_stress','threshold_stress_slip')
outputID = merge(thresholdstress_ID,undefined_ID,prm%sum_N_sl>0)
2018-11-28 10:29:03 +05:30
end select
2018-11-28 21:15:45 +05:30
2018-12-05 03:00:07 +05:30
if (outputID /= undefined_ID) then
2018-11-28 21:15:45 +05:30
plastic_disloUCLA_output(i,phase_plasticityInstance(p)) = outputs(i)
plastic_disloUCLA_sizePostResult(i,phase_plasticityInstance(p)) = prm%sum_N_sl
2018-11-28 21:15:45 +05:30
prm%outputID = [prm%outputID, outputID]
endif
2019-01-06 04:25:10 +05:30
enddo
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2019-01-07 12:34:02 +05:30
NipcMyPhase = count(material_phase == p)
sizeDotState = size(['rho_mob ','rho_dip ','gamma_sl']) * prm%sum_N_sl
sizeState = sizeDotState
2016-01-16 12:36:34 +05:30
2019-03-19 03:23:24 +05:30
call material_allocatePlasticState(p,NipcMyPhase,sizeState,sizeDotState,0, &
prm%sum_N_sl,0,0)
2018-12-05 03:00:07 +05:30
plasticState(p)%sizePostResults = sum(plastic_disloUCLA_sizePostResult(:,phase_plasticityInstance(p)))
!--------------------------------------------------------------------------------------------------
! locally defined state aliases and initialization of state0 and aTolState
2019-03-19 03:23:24 +05:30
startIndex = 1
endIndex = prm%sum_N_sl
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,:)
2018-12-05 03:00:07 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolRho
2016-04-25 02:06:35 +05:30
2019-03-19 03:23:24 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_sl
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,:)
2018-12-05 03:00:07 +05:30
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolRho
2016-04-25 02:06:35 +05:30
2019-03-19 03:23:24 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%sum_N_sl
2019-03-19 03:27:23 +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,:)
2016-04-25 02:06:35 +05:30
allocate(dst%mfp(prm%sum_N_sl,NipcMyPhase), source=0.0_pReal)
allocate(dst%dislocationSpacing(prm%sum_N_sl,NipcMyPhase),source=0.0_pReal)
allocate(dst%threshold_stress(prm%sum_N_sl,NipcMyPhase), source=0.0_pReal)
2018-12-21 20:31:16 +05:30
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
end associate
2018-12-04 04:36:46 +05:30
enddo
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
!--------------------------------------------------------------------------------------------------
2018-12-09 22:05:48 +05:30
pure subroutine plastic_disloUCLA_LpAndItsTangent(Lp,dLp_dMp,Mp,Temperature,instance,of)
2018-12-05 03:00:07 +05:30
implicit none
2018-12-21 20:31:16 +05:30
real(pReal), dimension(3,3), intent(out) :: &
Lp !< plastic velocity gradient
2018-12-21 20:31:16 +05:30
real(pReal), dimension(3,3,3,3), intent(out) :: &
dLp_dMp !< derivative of Lp with respect to the Mandel stress
real(pReal), dimension(3,3), intent(in) :: &
2018-12-21 20:31:16 +05:30
Mp !< Mandel stress
2019-03-19 03:27:23 +05:30
real(pReal), intent(in) :: &
2018-12-21 20:31:16 +05:30
temperature !< temperature
2019-03-19 03:27:23 +05:30
integer, intent(in) :: &
instance, &
of
2019-03-19 03:23:24 +05:30
integer :: &
i,k,l,m,n
real(pReal), dimension(param(instance)%sum_N_sl) :: &
gdot_pos,gdot_neg, &
dgdot_dtau_pos,dgdot_dtau_neg
2018-12-05 03:00:07 +05:30
Lp = 0.0_pReal
2018-11-28 11:14:32 +05:30
dLp_dMp = 0.0_pReal
2019-01-07 11:54:02 +05:30
associate(prm => param(instance))
2018-12-05 03:00:07 +05:30
call kinetics(Mp,Temperature,instance,of,gdot_pos,gdot_neg,dgdot_dtau_pos,dgdot_dtau_neg)
do i = 1, prm%sum_N_sl
Lp = Lp + (gdot_pos(i)+gdot_neg(i))*prm%Schmid(1:3,1:3,i)
2019-03-19 03:23:24 +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) &
+ dgdot_dtau_pos(i) * prm%Schmid(k,l,i) * prm%nonSchmid_pos(m,n,i) &
+ dgdot_dtau_neg(i) * prm%Schmid(k,l,i) * prm%nonSchmid_neg(m,n,i)
2019-01-07 12:34:02 +05:30
enddo
2019-01-07 11:54:02 +05:30
2018-12-05 02:03:32 +05:30
end associate
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief calculates the rate of change of microstructure
!--------------------------------------------------------------------------------------------------
2018-11-30 14:34:41 +05:30
subroutine plastic_disloUCLA_dotState(Mp,Temperature,instance,of)
use prec, only: &
2016-05-29 14:15:03 +05:30
tol_math_check, &
dEq0
use math, only: &
2018-12-05 01:20:02 +05:30
PI, &
math_clip
implicit none
2019-01-07 12:34:02 +05:30
real(pReal), dimension(3,3), intent(in) :: &
2018-12-21 17:16:43 +05:30
Mp !< Mandel stress
2019-01-07 12:34:02 +05:30
real(pReal), intent(in) :: &
2018-12-21 20:31:16 +05:30
temperature !< temperature
2019-03-19 03:23:24 +05:30
integer, intent(in) :: &
2019-01-07 12:34:02 +05:30
instance, &
of
real(pReal) :: &
2018-12-05 01:20:02 +05:30
VacancyDiffusion
real(pReal), dimension(param(instance)%sum_N_sl) :: &
gdot_pos, gdot_neg,&
tau_pos,&
tau_neg, &
2018-12-21 20:31:16 +05:30
DotRhoDipFormation, ClimbVelocity, EdgeDipDistance, &
2018-12-05 01:20:02 +05:30
DotRhoEdgeDipClimb
2016-04-25 02:06:35 +05:30
2018-12-05 03:00:07 +05:30
associate(prm => param(instance), stt => state(instance),dot => dotState(instance), dst => dependentState(instance))
2018-12-05 00:05:29 +05:30
2018-12-21 20:31:16 +05:30
call kinetics(Mp,Temperature,instance,of,&
gdot_pos,gdot_neg, &
2019-03-19 10:40:23 +05:30
tau_pos_out = tau_pos,tau_neg_out = tau_neg)
2018-12-05 03:00:07 +05:30
2019-03-19 03:27:23 +05:30
dot%gamma_sl(:,of) = (gdot_pos+gdot_neg) ! ToDo: needs to be abs
2018-12-05 02:03:32 +05:30
VacancyDiffusion = prm%D0*exp(-prm%Qsd/(kB*Temperature))
2018-12-05 01:20:02 +05:30
where(dEq0(tau_pos)) ! ToDo: use avg of pos and neg
2018-12-05 01:20:02 +05:30
DotRhoDipFormation = 0.0_pReal
DotRhoEdgeDipClimb = 0.0_pReal
else where
EdgeDipDistance = math_clip((3.0_pReal*prm%mu*prm%b_sl)/(16.0_pReal*PI*abs(tau_pos)), &
2018-12-21 20:31:16 +05:30
prm%minDipDistance, & ! lower limit
dst%mfp(:,of)) ! upper limit
DotRhoDipFormation = merge(((2.0_pReal*EdgeDipDistance)/prm%b_sl)* stt%rho_mob(:,of)*abs(dot%gamma_sl(:,of)), & ! ToDo: ignore region of spontaneous annihilation
2018-12-05 01:20:02 +05:30
0.0_pReal, &
prm%dipoleformation)
2018-12-09 19:19:08 +05:30
ClimbVelocity = (3.0_pReal*prm%mu*VacancyDiffusion*prm%atomicVolume/(2.0_pReal*pi*kB*Temperature)) &
2018-12-05 01:20:02 +05:30
* (1.0_pReal/(EdgeDipDistance+prm%minDipDistance))
DotRhoEdgeDipClimb = (4.0_pReal*ClimbVelocity*stt%rho_dip(:,of))/(EdgeDipDistance-prm%minDipDistance) ! ToDo: Discuss with Franz: Stress dependency?
2018-12-05 01:20:02 +05:30
end where
2018-11-28 00:19:04 +05:30
dot%rho_mob(:,of) = abs(dot%gamma_sl(:,of))/(prm%b_sl*dst%mfp(:,of)) & ! multiplication
2018-12-05 01:35:34 +05:30
- DotRhoDipFormation &
- (2.0_pReal*prm%minDipDistance)/prm%b_sl*stt%rho_mob(:,of)*abs(dot%gamma_sl(:,of)) !* Spontaneous annihilation of 2 single edge dislocations
dot%rho_dip(:,of) = DotRhoDipFormation &
- (2.0_pReal*prm%minDipDistance)/prm%b_sl*stt%rho_dip(:,of)*abs(dot%gamma_sl(:,of)) & !* Spontaneous annihilation of a single edge dislocation with a dipole constituent
2018-12-05 01:35:34 +05:30
- DotRhoEdgeDipClimb
2018-11-28 00:19:04 +05:30
2018-12-21 20:31:16 +05:30
end associate
2018-12-05 03:00:07 +05:30
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_dotState
2018-12-05 03:00:07 +05:30
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine plastic_disloUCLA_dependentState(instance,of)
implicit none
2019-03-19 03:23:24 +05:30
integer, intent(in) :: &
2019-01-07 12:34:02 +05:30
instance, &
of
2019-01-07 12:06:11 +05:30
2019-03-19 03:23:24 +05:30
integer :: &
2019-01-07 12:06:11 +05:30
i
associate(prm => param(instance), stt => state(instance),dst => dependentState(instance))
forall (i = 1:prm%sum_N_sl)
dst%dislocationSpacing(i,of) = sqrt(dot_product(stt%rho_mob(:,of)+stt%rho_dip(:,of), &
2019-01-07 12:06:11 +05:30
prm%forestProjectionEdge(:,i)))
dst%threshold_stress(i,of) = prm%mu*prm%b_sl(i) &
* sqrt(dot_product(stt%rho_mob(:,of)+stt%rho_dip(:,of), &
2019-03-19 10:40:23 +05:30
prm%h_sl_sl(:,i)))
2019-01-07 12:06:11 +05:30
end forall
dst%mfp(:,of) = prm%D/(1.0_pReal+prm%D*dst%dislocationSpacing(:,of)/prm%Clambda)
2019-01-07 12:06:11 +05:30
dst%dislocationSpacing(:,of) = dst%mfp(:,of) ! ToDo: Hack to recover wrong behavior for the moment
end associate
end subroutine plastic_disloUCLA_dependentState
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2018-11-30 14:34:41 +05:30
function plastic_disloUCLA_postResults(Mp,Temperature,instance,of) result(postResults)
use prec, only: &
2017-11-21 04:13:06 +05:30
dEq, dNeq0
use math, only: &
PI, &
2018-12-03 15:18:37 +05:30
math_mul33xx33
implicit none
2018-12-03 15:18:37 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
real(pReal), intent(in) :: &
2018-12-21 20:31:16 +05:30
temperature !< temperature
2019-03-19 03:23:24 +05:30
integer, intent(in) :: &
2018-12-03 15:18:37 +05:30
instance, &
of
2018-12-03 15:18:37 +05:30
real(pReal), dimension(sum(plastic_disloUCLA_sizePostResult(:,instance))) :: &
postResults
2019-03-19 03:23:24 +05:30
integer :: &
2018-12-03 15:18:37 +05:30
o,c,i
real(pReal), dimension(param(instance)%sum_N_sl) :: &
gdot_pos,gdot_neg
2018-11-30 14:34:41 +05:30
2019-03-19 03:23:24 +05:30
c = 0
2018-12-03 15:18:37 +05:30
associate(prm => param(instance), stt => state(instance), dst => dependentState(instance))
2019-03-19 03:23:24 +05:30
outputsLoop: do o = 1,size(prm%outputID)
2018-12-03 15:18:37 +05:30
select case(prm%outputID(o))
2019-03-19 10:16:47 +05:30
case (rho_mob_ID)
postResults(c+1:c+prm%sum_N_sl) = stt%rho_mob(1:prm%sum_N_sl,of)
2019-03-19 10:16:47 +05:30
case (rho_dip_ID)
postResults(c+1:c+prm%sum_N_sl) = stt%rho_dip(1:prm%sum_N_sl,of)
2019-03-19 10:16:47 +05:30
case (gamma_dot_sl_ID)
call kinetics(Mp,Temperature,instance,of,gdot_pos,gdot_neg)
postResults(c+1:c+prm%sum_N_sl) = gdot_pos + gdot_neg
2019-03-19 10:16:47 +05:30
case (gamma_sl_ID)
postResults(c+1:c+prm%sum_N_sl) = stt%gamma_sl(1:prm%sum_N_sl, of)
2019-03-19 10:40:23 +05:30
case (Lambda_sl_ID)
postResults(c+1:c+prm%sum_N_sl) = dst%mfp(1:prm%sum_N_sl, of)
2018-12-03 15:18:37 +05:30
case (thresholdstress_ID)
postResults(c+1:c+prm%sum_N_sl) = dst%threshold_stress(1:prm%sum_N_sl,of)
2018-12-03 15:18:37 +05:30
end select
c = c + prm%sum_N_sl
2018-12-03 15:18:37 +05:30
enddo outputsLoop
2018-12-03 15:18:37 +05:30
end associate
end function plastic_disloUCLA_postResults
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine plastic_disloUCLA_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 03:23:24 +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_disloUCLA_results
2019-01-07 11:54:02 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Shear rates on slip systems, their derivatives with respect to resolved stress and the
! resolved stresss
2019-01-07 11:54:02 +05:30
!> @details Derivatives and resolved stress are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end
!--------------------------------------------------------------------------------------------------
2018-12-21 17:32:51 +05:30
pure subroutine kinetics(Mp,Temperature,instance,of, &
2019-03-19 10:40:23 +05:30
gamma_pos,gamma_neg,dgamma_dtau_pos,dgamma_dtau_neg,tau_pos_out,tau_neg_out)
use prec, only: &
tol_math_check, &
dEq, dNeq0
use math, only: &
2018-12-21 17:16:43 +05:30
PI, &
math_mul33xx33
implicit none
2018-11-28 11:14:32 +05:30
real(pReal), dimension(3,3), intent(in) :: &
2018-12-21 20:31:16 +05:30
Mp !< Mandel stress
real(pReal), intent(in) :: &
2018-12-21 20:31:16 +05:30
temperature !< temperature
2019-03-19 03:23:24 +05:30
integer, intent(in) :: &
instance, &
of
real(pReal), intent(out), dimension(param(instance)%sum_N_sl) :: &
gamma_pos, &
gamma_neg
real(pReal), intent(out), optional, dimension(param(instance)%sum_N_sl) :: &
dgamma_dtau_pos, &
dgamma_dtau_neg, &
2019-03-19 10:40:23 +05:30
tau_pos_out, &
tau_neg_out
real(pReal), dimension(param(instance)%sum_N_sl) :: &
2018-12-21 20:31:16 +05:30
StressRatio, &
2018-12-09 19:55:54 +05:30
StressRatio_p,StressRatio_pminus1, &
dvel, vel, &
tau_pos,tau_neg, &
t_n, t_k, dtk,dtn, &
2018-12-21 20:31:16 +05:30
needsGoodName ! ToDo: @Karo: any idea?
2019-03-19 03:23:24 +05:30
integer :: j
2018-12-21 20:31:16 +05:30
associate(prm => param(instance), stt => state(instance), dst => dependentState(instance))
2019-01-07 11:54:02 +05:30
do j = 1, prm%sum_N_sl
tau_pos(j) = math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,j))
tau_neg(j) = math_mul33xx33(Mp,prm%nonSchmid_neg(1:3,1:3,j))
2018-12-09 19:19:08 +05:30
enddo
2019-01-07 11:54:02 +05:30
2019-03-19 10:40:23 +05:30
if (present(tau_pos_out)) tau_pos_out = tau_pos
if (present(tau_neg_out)) tau_neg_out = tau_neg
2018-12-05 03:00:07 +05:30
2019-03-19 10:40:23 +05:30
associate(BoltzmannRatio => prm%delta_F/(kB*Temperature), &
DotGamma0 => stt%rho_mob(:,of)*prm%b_sl*prm%v0, &
2018-12-21 21:34:26 +05:30
effectiveLength => dst%mfp(:,of) - prm%w)
2018-12-09 19:55:54 +05:30
significantPositiveTau: where(abs(tau_pos)-dst%threshold_stress(:,of) > tol_math_check)
2019-03-19 10:40:23 +05:30
StressRatio = (abs(tau_pos)-dst%threshold_stress(:,of))/prm%tau_0
2018-12-21 17:16:43 +05:30
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
2018-12-21 20:31:16 +05:30
needsGoodName = exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)
t_n = prm%b_sl/(needsGoodName*prm%omega*effectiveLength)
t_k = effectiveLength * prm%B /(2.0_pReal*prm%b_sl*tau_pos) ! our definition of tk is different with the one in dislotwin
2018-12-21 17:16:43 +05:30
vel = prm%kink_height/(t_n + t_k)
2018-12-05 02:03:32 +05:30
gamma_pos = DotGamma0 * sign(vel,tau_pos) * 0.5_pReal
2018-12-21 17:32:51 +05:30
else where significantPositiveTau
gamma_pos = 0.0_pReal
2018-12-21 17:32:51 +05:30
end where significantPositiveTau
2018-12-21 17:16:43 +05:30
if (present(dgamma_dtau_pos)) then
significantPositiveTau2: where(abs(tau_pos)-dst%threshold_stress(:,of) > tol_math_check)
dtn = t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
2019-03-19 10:40:23 +05:30
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_0
dtk = t_k / tau_pos
dvel = prm%kink_height * (dtk + dtn) / (t_n + t_k)**2.0_pReal
dgamma_dtau_pos = DotGamma0 * dvel* 0.5_pReal
2018-12-21 17:32:51 +05:30
else where significantPositiveTau2
dgamma_dtau_pos = 0.0_pReal
2018-12-21 17:32:51 +05:30
end where significantPositiveTau2
2018-12-21 20:31:16 +05:30
endif
2018-12-21 17:16:43 +05:30
significantNegativeTau: where(abs(tau_neg)-dst%threshold_stress(:,of) > tol_math_check)
2019-03-19 10:40:23 +05:30
StressRatio = (abs(tau_neg)-dst%threshold_stress(:,of))/prm%tau_0
2018-12-21 17:16:43 +05:30
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
2018-12-21 20:31:16 +05:30
needsGoodName = exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)
2018-12-21 17:16:43 +05:30
t_n = prm%b_sl/(needsGoodName*prm%omega*effectiveLength)
t_k = effectiveLength * prm%B /(2.0_pReal*prm%b_sl*tau_pos) ! our definition of tk is different with the one in dislotwin
vel = prm%kink_height/(t_n + t_k)
2018-12-05 02:03:32 +05:30
gamma_neg = DotGamma0 * sign(vel,tau_neg) * 0.5_pReal
2018-12-21 17:32:51 +05:30
else where significantNegativeTau
gamma_neg = 0.0_pReal
2018-12-21 17:32:51 +05:30
end where significantNegativeTau
2018-12-21 17:16:43 +05:30
if (present(dgamma_dtau_neg)) then
significantNegativeTau2: where(abs(tau_neg)-dst%threshold_stress(:,of) > tol_math_check)
dtn = t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
2019-03-19 10:40:23 +05:30
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_0
dtk = t_k / tau_pos
dvel = prm%kink_height * (dtk + dtn) / (t_n + t_k)**2.0_pReal
dgamma_dtau_neg = DotGamma0 * dvel * 0.5_pReal
2018-12-21 17:32:51 +05:30
else where significantNegativeTau2
dgamma_dtau_neg = 0.0_pReal
2018-12-21 17:32:51 +05:30
end where significantNegativeTau2
2018-12-21 20:31:16 +05:30
end if
2019-03-19 10:40:23 +05:30
2018-12-21 20:31:16 +05:30
end associate
2018-12-21 17:32:51 +05:30
end associate
2018-12-05 02:03:32 +05:30
end subroutine kinetics
2015-01-21 20:44:00 +05:30
end module plastic_disloUCLA