DAMASK_EICMD/src/plastic_disloUCLA.f90

927 lines
43 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 material subroutine incoprorating dislocation and twinning physics
!> @details to be done
!--------------------------------------------------------------------------------------------------
2015-01-21 20:44:00 +05:30
module plastic_disloUCLA
use prec, only: &
pReal, &
pInt
implicit none
private
integer(pInt), dimension(:,:), allocatable, target, public :: &
2016-01-16 12:36:34 +05:30
plastic_disloUCLA_sizePostResult !< size of each post result output
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
2016-01-16 12:36:34 +05:30
integer(pInt), dimension(:), allocatable, private :: &
2016-04-25 02:06:35 +05:30
plastic_disloUCLA_totalNslip !< total number of active slip systems for each instance
integer(pInt), dimension(:,:), allocatable, private :: &
2016-04-25 02:06:35 +05:30
plastic_disloUCLA_Nslip !< number of active slip systems for each family and instance
real(pReal), dimension(:), allocatable, private :: &
2016-01-16 12:36:34 +05:30
plastic_disloUCLA_CAtomicVolume, & !< atomic volume in Bugers vector unit
plastic_disloUCLA_Qsd, & !< activation energy for dislocation climb
plastic_disloUCLA_CEdgeDipMinDistance, & !<
plastic_disloUCLA_dipoleFormationFactor !< scaling factor for dipole formation: 0: off, 1: on. other values not useful
real(pReal), dimension(:,:,:), allocatable, private :: &
2016-01-16 12:36:34 +05:30
plastic_disloUCLA_forestProjectionEdge !< matrix of forest projections of edge dislocations for each instance
enum, bind(c)
enumerator :: undefined_ID, &
2018-11-28 10:29:03 +05:30
rho_ID, &
rhoDip_ID, &
shearrate_ID, &
accumulatedshear_ID, &
mfp_ID, &
resolvedstress_ID, &
thresholdstress_ID, &
dipoledistance_ID, &
stressexponent_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, &
2018-11-29 15:01:02 +05:30
grainSize, &
2018-11-30 12:55:23 +05:30
SolidSolutionStrength, & !< Strength due to elements in solid solution
2018-12-03 15:55:29 +05:30
mu, &
D0 !< prefactor for self-diffusion coefficient
2018-11-28 00:30:45 +05:30
real(pReal), allocatable, dimension(:) :: &
2018-11-29 15:07:06 +05:30
B, & !< friction coeff. B (kMC)
2018-11-28 10:29:03 +05:30
rho0, & !< initial edge dislocation density per slip system for each family and instance
rhoDip0, & !< initial edge dipole density per slip system for each family and instance
burgers, & !< absolute length of burgers vector [m] for each slip system and instance
H0kp, & !< activation energy for glide [J] for each slip system and instance
v0, & !< dislocation velocity prefactor [m/s] for each family and instance
CLambda, & !< Adj. parameter for distance between 2 forest dislocations for each slip system and instance
p, & !< p-exponent in glide velocity
q, & !< q-exponent in glide velocity
!* mobility law parameters
kink_height, & !< height of the kink pair
kink_width, & !< width of the kink pair
omega, & !< attempt frequency for kink pair nucleation
2018-11-28 10:29:03 +05:30
viscosity, & !< friction coeff. B (kMC)
!*
2018-11-29 13:14:31 +05:30
tau_Peierls, &
2018-11-28 00:30:45 +05:30
nonSchmidCoeff
real(pReal), allocatable, dimension(:,:) :: &
interaction_SlipSlip !< slip resistance from slip activity
real(pReal), allocatable, dimension(:,:,:) :: &
Schmid_slip, &
Schmid_twin, &
nonSchmid_pos, &
nonSchmid_neg
integer(pInt) :: &
totalNslip !< total number of active slip system
integer(pInt), allocatable, dimension(:) :: &
Nslip !< number of active slip systems for each family
integer(kind(undefined_ID)), allocatable, dimension(:) :: &
outputID !< ID of each post result output
end type !< container type for internal constitutive parameters
type(tParameters), dimension(:), allocatable, private :: param !< containers of constitutive parameters (len Ninstance)
type, private :: tDisloUCLAState
real(pReal), pointer, dimension(:,:) :: &
rhoEdge, &
rhoEdgeDip, &
2018-11-30 14:34:41 +05:30
accshear_slip, &
whole
end type
type, private :: tDisloUCLAMicrostructure
real(pReal), allocatable, dimension(:,:) :: &
mfp, &
threshold_stress
end type tDisloUCLAMicrostructure
type(tDisloUCLAState ), allocatable, dimension(:), private :: &
state, &
dotState
type(tDisloUCLAMicrostructure), allocatable, dimension(:), private :: &
microstructure
public :: &
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_init, &
plastic_disloUCLA_microstructure, &
plastic_disloUCLA_LpAndItsTangent, &
plastic_disloUCLA_dotState, &
plastic_disloUCLA_postResults
private :: &
2018-11-29 13:51:58 +05:30
kinetics
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2015-01-21 20:44:00 +05:30
subroutine plastic_disloUCLA_init(fileUnit)
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
use, intrinsic :: iso_fortran_env, only: &
compiler_version, &
compiler_options
#endif
use debug, only: &
debug_level,&
debug_constitutive,&
debug_levelBasic
use math, only: &
math_Mandel3333to66, &
math_Voigt66to3333, &
math_mul3x3, &
math_expand
use IO, only: &
IO_error, &
2018-11-30 13:06:56 +05:30
IO_timeStamp
use material, only: &
phase_plasticity, &
phase_plasticityInstance, &
phase_Noutput, &
2015-01-15 16:26:15 +05:30
PLASTICITY_DISLOUCLA_label, &
PLASTICITY_DISLOUCLA_ID, &
material_phase, &
2018-11-27 23:58:00 +05:30
plasticState, &
material_allocatePlasticState
use config, only: &
2018-11-28 00:30:45 +05:30
MATERIAL_partPhase, &
config_phase
use lattice
implicit none
integer(pInt), intent(in) :: fileUnit
2018-11-29 13:51:58 +05:30
integer(pInt) :: maxNinstance,phase,maxTotalNslip,&
2018-11-28 10:29:03 +05:30
f,instance,j,k,o,ns, i, &
2018-11-30 13:06:56 +05:30
outputSize, &
offset_slip, index_myFamily, index_otherFamily, &
2018-11-28 00:30:45 +05:30
startIndex, endIndex, p
2018-11-29 15:01:02 +05:30
integer(pInt) :: sizeState, sizeDotState
integer(pInt) :: NofMyPhase
character(len=65536) :: &
2018-12-03 15:55:29 +05:30
structure = ''
2018-11-28 10:29:03 +05:30
character(len=65536), dimension(:), allocatable :: outputs
integer(kind(undefined_ID)) :: outputID
2018-11-28 00:30:45 +05:30
integer(pInt), dimension(0), parameter :: emptyIntArray = [integer(pInt)::]
real(pReal), dimension(0), parameter :: emptyRealArray = [real(pReal)::]
character(len=65536), dimension(0), parameter :: emptyStringArray = [character(len=65536)::]
write(6,'(/,a)') ' <<<+- constitutive_'//PLASTICITY_DISLOUCLA_label//' init -+>>>'
write(6,'(/,a)') ' Cereceda et al., International Journal of Plasticity 78, 2016, 242-256'
2018-04-09 18:34:37 +05:30
write(6,'(/,a)') ' http://dx.doi.org/10.1016/j.ijplas.2015.09.002'
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
#include "compilation_info.f90"
2015-01-15 16:26:15 +05:30
maxNinstance = int(count(phase_plasticity == PLASTICITY_DISLOUCLA_ID),pInt)
if (maxNinstance == 0_pInt) return
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0_pInt) &
write(6,'(a16,1x,i5,/)') '# instances:',maxNinstance
2015-01-21 20:44:00 +05:30
allocate(plastic_disloUCLA_sizePostResult(maxval(phase_Noutput),maxNinstance),source=0_pInt)
allocate(plastic_disloUCLA_output(maxval(phase_Noutput),maxNinstance))
plastic_disloUCLA_output = ''
2018-12-03 15:55:29 +05:30
2015-01-21 20:44:00 +05:30
allocate(plastic_disloUCLA_Nslip(lattice_maxNslipFamily,maxNinstance), source=0_pInt)
allocate(plastic_disloUCLA_totalNslip(maxNinstance), source=0_pInt)
allocate(plastic_disloUCLA_CAtomicVolume(maxNinstance), source=0.0_pReal)
2018-12-03 15:55:29 +05:30
2015-01-21 20:44:00 +05:30
allocate(plastic_disloUCLA_Qsd(maxNinstance), source=0.0_pReal)
allocate(plastic_disloUCLA_CEdgeDipMinDistance(maxNinstance), source=0.0_pReal)
allocate(plastic_disloUCLA_dipoleFormationFactor(maxNinstance), source=1.0_pReal) !should be on by default
2018-11-28 00:30:45 +05:30
allocate(param(maxNinstance))
allocate(state(maxNinstance))
allocate(dotState(maxNinstance))
allocate(microstructure(maxNinstance))
2018-11-28 00:30:45 +05:30
do p = 1_pInt, size(phase_plasticityInstance)
if (phase_plasticity(p) /= PLASTICITY_DISLOUCLA_ID) cycle
associate(prm => param(phase_plasticityInstance(p)), &
dot => dotState(phase_plasticityInstance(p)), &
stt => state(phase_plasticityInstance(p)))
structure = config_phase(p)%getString('lattice_structure')
2018-11-30 12:55:23 +05:30
prm%mu = lattice_mu(p)
2018-11-28 00:30:45 +05:30
prm%aTolRho = config_phase(p)%getFloat('atol_rho')
2018-11-28 00:30:45 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
prm%Nslip = config_phase(p)%getInts('nslip',defaultVal=emptyIntArray)
prm%totalNslip = sum(prm%Nslip)
slipActive: if (prm%totalNslip > 0_pInt) then
prm%Schmid_slip = lattice_SchmidMatrix_slip(prm%Nslip,structure(1:3),&
config_phase(p)%getFloat('c/a',defaultVal=0.0_pReal))
if(structure=='bcc') then
prm%nonSchmidCoeff = config_phase(p)%getFloats('nonschmid_coefficients',&
defaultVal = emptyRealArray)
prm%nonSchmid_pos = lattice_nonSchmidMatrix(prm%Nslip,prm%nonSchmidCoeff,+1_pInt)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(prm%Nslip,prm%nonSchmidCoeff,-1_pInt)
else
prm%nonSchmid_pos = prm%Schmid_slip
prm%nonSchmid_neg = prm%Schmid_slip
endif
prm%interaction_SlipSlip = lattice_interaction_SlipSlip(prm%Nslip, &
config_phase(p)%getFloats('interaction_slipslip'), &
structure(1:3))
prm%rho0 = config_phase(p)%getFloats('rhoedge0')
prm%rhoDip0 = config_phase(p)%getFloats('rhoedgedip0')
prm%burgers = config_phase(p)%getFloats('slipburgers')
prm%H0kp = config_phase(p)%getFloats('qedge')
2018-11-29 15:01:02 +05:30
prm%v0 = config_phase(p)%getFloats('v0')
prm%clambda = config_phase(p)%getFloats('clambdaslip')
2018-11-29 13:14:31 +05:30
prm%tau_Peierls = config_phase(p)%getFloats('tau_peierls')
prm%p = config_phase(p)%getFloats('p_slip',defaultVal=[(1.0_pReal,i=1_pInt,size(prm%Nslip))])
prm%q = config_phase(p)%getFloats('q_slip',defaultVal=[(1.0_pReal,i=1_pInt,size(prm%Nslip))])
prm%kink_height = config_phase(p)%getFloats('kink_height')
prm%kink_width = config_phase(p)%getFloats('kink_width')
prm%omega = config_phase(p)%getFloats('omega')
2018-11-29 15:07:06 +05:30
prm%B = config_phase(p)%getFloats('friction_coeff')
2018-11-28 10:29:03 +05:30
!prm%viscosity = config_phase(p)%getFloats('viscosity')
2018-11-29 15:01:02 +05:30
prm%SolidSolutionStrength = config_phase(p)%getFloat('solidsolutionstrength')
2018-11-29 13:14:31 +05:30
prm%grainSize = config_phase(p)%getFloat('grainsize')
2018-12-03 15:55:29 +05:30
prm%D0 = config_phase(p)%getFloat('d0')
2018-11-29 15:01:02 +05:30
plastic_disloUCLA_Qsd(phase_plasticityInstance(p)) = config_phase(p)%getFloat('qsd')
plastic_disloUCLA_CEdgeDipMinDistance(phase_plasticityInstance(p)) = config_phase(p)%getFloat('cedgedipmindistance')
plastic_disloUCLA_CAtomicVolume(phase_plasticityInstance(p)) = config_phase(p)%getFloat('catomicvolume')
plastic_disloUCLA_dipoleFormationFactor(phase_plasticityInstance(p)) = config_phase(p)%getFloat('dipoleformationfactor')
2018-11-29 13:14:31 +05:30
! expand: family => system
prm%rho0 = math_expand(prm%rho0, prm%Nslip)
prm%rhoDip0 = math_expand(prm%rhoDip0, prm%Nslip)
prm%q = math_expand(prm%q, prm%Nslip)
prm%p = math_expand(prm%p, prm%Nslip)
prm%H0kp = math_expand(prm%H0kp, prm%Nslip)
prm%burgers = math_expand(prm%burgers, prm%Nslip)
prm%kink_height = math_expand(prm%kink_height, prm%Nslip)
prm%kink_width = math_expand(prm%kink_width, prm%Nslip)
prm%omega = math_expand(prm%omega, prm%Nslip)
2018-11-29 13:14:31 +05:30
prm%tau_Peierls = math_expand(prm%tau_Peierls, prm%Nslip)
2018-11-29 15:01:02 +05:30
prm%v0 = math_expand(prm%v0, prm%Nslip)
2018-11-29 15:07:06 +05:30
prm%B = math_expand(prm%B, prm%Nslip)
prm%clambda = math_expand(prm%clambda, prm%Nslip)
instance = phase_plasticityInstance(p)
2018-11-30 12:55:23 +05:30
plastic_disloUCLA_totalNslip(instance) = prm%totalNslip
2018-12-03 15:55:29 +05:30
!if (plastic_disloUCLA_CAtomicVolume(instance) <= 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='cAtomicVolume ('//PLASTICITY_DISLOUCLA_label//')')
if (prm%D0 <= 0.0_pReal) &
call IO_error(211_pInt,el=instance,ext_msg='D0 ('//PLASTICITY_DISLOUCLA_label//')')
if (plastic_disloUCLA_Qsd(instance) <= 0.0_pReal) &
call IO_error(211_pInt,el=instance,ext_msg='Qsd ('//PLASTICITY_DISLOUCLA_label//')')
! if (plastic_disloUCLA_aTolRho(instance) <= 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='aTolRho ('//PLASTICITY_DISLOUCLA_label//')')
else slipActive
allocate(prm%rho0(0))
allocate(prm%rhoDip0(0))
2018-11-28 00:30:45 +05:30
endif slipActive
2018-11-28 10:29:03 +05:30
!--------------------------------------------------------------------------------------------------
! phase outputs
#if defined(__GFORTRAN__)
outputs = ['GfortranBug86277']
outputs = config_phase(p)%getStrings('(output)',defaultVal=outputs)
if (outputs(1) == 'GfortranBug86277') outputs = emptyStringArray
#else
outputs = config_phase(p)%getStrings('(output)',defaultVal=emptyStringArray)
#endif
allocate(prm%outputID(0))
do i = 1_pInt, size(outputs)
outputID = undefined_ID
outputSize = prm%totalNslip
select case(trim(outputs(i)))
case ('edge_density')
2018-11-28 21:15:45 +05:30
outputID = merge(rho_ID,undefined_ID,prm%totalNslip>0_pInt)
2018-11-28 10:29:03 +05:30
case ('dipole_density')
2018-11-28 21:15:45 +05:30
outputID = merge(rhoDip_ID,undefined_ID,prm%totalNslip>0_pInt)
case ('shear_rate','shearrate','shear_rate_slip','shearrate_slip')
outputID = merge(shearrate_ID,undefined_ID,prm%totalNslip>0_pInt)
case ('accumulated_shear','accumulatedshear','accumulated_shear_slip')
outputID = merge(accumulatedshear_ID,undefined_ID,prm%totalNslip>0_pInt)
case ('mfp','mfp_slip')
outputID = merge(mfp_ID,undefined_ID,prm%totalNslip>0_pInt)
case ('resolved_stress','resolved_stress_slip')
outputID = merge(resolvedstress_ID,undefined_ID,prm%totalNslip>0_pInt)
case ('threshold_stress','threshold_stress_slip')
outputID = merge(thresholdstress_ID,undefined_ID,prm%totalNslip>0_pInt)
2018-11-28 10:29:03 +05:30
case ('edge_dipole_distance')
2018-11-28 21:15:45 +05:30
outputID = merge(dipoleDistance_ID,undefined_ID,prm%totalNslip>0_pInt)
2018-11-28 10:29:03 +05:30
case ('stress_exponent')
2018-11-28 21:15:45 +05:30
outputID = merge(stressexponent_ID,undefined_ID,prm%totalNslip>0_pInt)
2018-11-28 10:29:03 +05:30
end select
2018-11-28 21:15:45 +05:30
if (outputID /= undefined_ID) then
plastic_disloUCLA_output(i,phase_plasticityInstance(p)) = outputs(i)
plastic_disloUCLA_sizePostResult(i,phase_plasticityInstance(p)) = outputSize
prm%outputID = [prm%outputID, outputID]
endif
2018-11-28 10:29:03 +05:30
enddo
2018-11-28 00:30:45 +05:30
end associate
enddo
2018-11-30 13:06:56 +05:30
!if (plastic_disloUCLA_rhoEdge0(f,instance) < 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='rhoEdge0 ('//PLASTICITY_DISLOUCLA_label//')')
!if (plastic_disloUCLA_rhoEdgeDip0(f,instance) < 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='rhoEdgeDip0 ('//PLASTICITY_DISLOUCLA_label//')')
!if (plastic_disloUCLA_burgersPerSlipFamily(f,instance) <= 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='slipBurgers ('//PLASTICITY_DISLOUCLA_label//')')
2018-11-29 15:01:02 +05:30
!if (plastic_disloUCLA_v0PerSlipFamily(f,instance) <= 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='v0 ('//PLASTICITY_DISLOUCLA_label//')')
2018-11-29 13:14:31 +05:30
!if (plastic_disloUCLA_tau_peierlsPerSlipFamily(f,instance) < 0.0_pReal) &
! call IO_error(211_pInt,el=instance,ext_msg='tau_peierls ('//PLASTICITY_DISLOUCLA_label//')')
2018-11-30 13:06:56 +05:30
!--------------------------------------------------------------------------------------------------
! allocation of variables whose size depends on the total number of active slip systems
2015-01-21 20:44:00 +05:30
maxTotalNslip = maxval(plastic_disloUCLA_totalNslip)
allocate(plastic_disloUCLA_forestProjectionEdge(maxTotalNslip,maxTotalNslip,maxNinstance), &
source=0.0_pReal)
2018-11-28 00:30:45 +05:30
initializeInstances: do phase = 1_pInt, size(phase_plasticity)
2016-04-25 02:06:35 +05:30
myPhase2: if (phase_plasticity(phase) == PLASTICITY_disloUCLA_ID) then
p = phase
NofMyPhase=count(material_phase==phase)
instance = phase_plasticityInstance(phase)
2015-01-21 20:44:00 +05:30
ns = plastic_disloUCLA_totalNslip(instance)
associate(prm => param(instance), stt=>state(instance),mse => microstructure(phase_plasticityInstance(p)))
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2016-01-16 12:36:34 +05:30
2016-04-25 02:06:35 +05:30
sizeDotState = int(size(['rhoEdge ','rhoEdgeDip ','accshearslip']),pInt) * ns
sizeState = sizeDotState
2016-01-16 12:36:34 +05:30
2018-11-27 23:58:00 +05:30
call material_allocatePlasticState(phase,NofMyPhase,sizeState,sizeDotState,0_pInt, &
ns,0_pInt,0_pInt)
2018-12-03 15:55:29 +05:30
plasticState(phase)%sizePostResults = sum(plastic_disloUCLA_sizePostResult(:,phase_plasticityInstance(p)))
2015-01-23 18:38:25 +05:30
offset_slip = 2_pInt*plasticState(phase)%nSlip
plasticState(phase)%slipRate => &
plasticState(phase)%dotState(offset_slip+1:offset_slip+plasticState(phase)%nSlip,1:NofMyPhase)
plasticState(phase)%accumulatedSlip => &
plasticState(phase)%state (offset_slip+1:offset_slip+plasticState(phase)%nSlip,1:NofMyPhase)
!* Process slip related parameters ------------------------------------------------
2018-11-29 15:07:06 +05:30
mySlipFamilies: do f = 1_pInt,size(prm%Nslip,1)
2018-11-30 12:55:23 +05:30
index_myFamily = sum(prm%Nslip(1:f-1_pInt)) ! index in truncated slip system list
mySlipSystems: do j = 1_pInt,prm%Nslip(f)
!* Calculation of forest projections for edge dislocations
2018-11-29 15:07:06 +05:30
otherSlipFamilies: do o = 1_pInt,size(prm%Nslip,1)
2018-11-30 12:55:23 +05:30
index_otherFamily = sum(prm%Nslip(1:o-1_pInt))
otherSlipSystems: do k = 1_pInt,prm%Nslip(o)
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_forestProjectionEdge(index_myFamily+j,index_otherFamily+k,instance) = &
abs(math_mul3x3(lattice_sn(:,sum(lattice_NslipSystem(1:f-1,phase))+j,phase), &
lattice_st(:,sum(lattice_NslipSystem(1:o-1,phase))+k,phase)))
enddo otherSlipSystems; enddo otherSlipFamilies
2016-04-25 02:06:35 +05:30
enddo mySlipSystems
enddo mySlipFamilies
2016-04-25 02:06:35 +05:30
startIndex=1_pInt
endIndex=ns
2018-11-29 11:57:35 +05:30
stt%rhoEdge=>plasticState(phase)%state(startIndex:endIndex,:)
stt%rhoEdge= spread(prm%rho0,2,NofMyPhase)
dotState(instance)%rhoEdge=>plasticState(phase)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolRho
2016-04-25 02:06:35 +05:30
startIndex=endIndex+1_pInt
endIndex=endIndex+ns
2018-11-29 11:57:35 +05:30
stt%rhoEdgeDip=>plasticState(phase)%state(startIndex:endIndex,:)
stt%rhoEdgeDip= spread(prm%rhoDip0,2,NofMyPhase)
2016-04-25 02:06:35 +05:30
dotState(instance)%rhoEdgeDip=>plasticState(phase)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolRho
2016-04-25 02:06:35 +05:30
startIndex=endIndex+1_pInt
endIndex=endIndex+ns
2018-11-29 11:57:35 +05:30
stt%accshear_slip=>plasticState(phase)%state(startIndex:endIndex,:)
2016-04-25 02:06:35 +05:30
dotState(instance)%accshear_slip=>plasticState(phase)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = 1e6_pReal
2016-04-25 02:06:35 +05:30
2018-11-30 14:34:41 +05:30
dotState(instance)%whole => plasticState(phase)%dotState
allocate(mse%mfp(prm%totalNslip,NofMyPhase),source=0.0_pReal)
allocate(mse%threshold_stress(prm%totalNslip,NofMyPhase),source=0.0_pReal)
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
end associate
endif myPhase2
enddo initializeInstances
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
2015-01-21 20:44:00 +05:30
subroutine plastic_disloUCLA_microstructure(temperature,ipc,ip,el)
use math, only: &
pi
use material, only: &
phase_plasticityInstance, &
phaseAt, phasememberAt, &
material_phase
implicit none
integer(pInt), intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
real(pReal), intent(in) :: &
temperature !< temperature at IP
integer(pInt) :: &
instance, &
2016-04-25 02:06:35 +05:30
ns,s, &
of
real(pReal), dimension(plastic_disloUCLA_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))) :: &
invLambdaSlip
!* Shortened notation
of = phasememberAt(ipc,ip,el)
2018-11-30 12:16:26 +05:30
instance = phase_plasticityInstance(phaseAt(ipc,ip,el))
2015-01-21 20:44:00 +05:30
ns = plastic_disloUCLA_totalNslip(instance)
2018-11-29 13:14:31 +05:30
associate(prm => param(instance), stt => state(instance),mse => microstructure(instance))
!* 1/mean free distance between 2 forest dislocations seen by a moving dislocation
forall (s = 1_pInt:ns) &
invLambdaSlip(s) = &
2018-11-29 11:57:35 +05:30
sqrt(dot_product((stt%rhoEdge(1_pInt:ns,of)+stt%rhoEdgeDip(1_pInt:ns,of)),&
2015-01-21 20:44:00 +05:30
plastic_disloUCLA_forestProjectionEdge(1:ns,s,instance)))/ &
prm%Clambda(s)
!* mean free path between 2 obstacles seen by a moving dislocation
2018-11-29 15:01:02 +05:30
mse%mfp(:,of) = prm%grainSize/(1.0_pReal+prm%grainSize*invLambdaSlip)
!* threshold stress for dislocation motion
forall (s = 1_pInt:ns) &
mse%threshold_stress(s,of) = &
2018-11-30 12:55:23 +05:30
prm%mu*prm%burgers(s)*&
2018-11-29 15:01:02 +05:30
sqrt(dot_product(stt%rhoEdge(1_pInt:ns,of)+stt%rhoEdgeDip(1_pInt:ns,of),&
prm%interaction_SlipSlip(s,1:ns)))
end associate
2018-11-29 13:14:31 +05:30
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_microstructure
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
!--------------------------------------------------------------------------------------------------
2018-11-28 11:14:32 +05:30
subroutine plastic_disloUCLA_LpAndItsTangent(Lp,dLp_dMp,Mp,Temperature,ipc,ip,el)
use material, only: &
material_phase, &
phase_plasticityInstance, &
phaseAt, phasememberAt
implicit none
integer(pInt), intent(in) :: ipc,ip,el
real(pReal), intent(in) :: Temperature
2018-11-28 11:14:32 +05:30
real(pReal), dimension(3,3), intent(in) :: Mp
real(pReal), dimension(3,3), intent(out) :: Lp
2018-11-28 11:14:32 +05:30
real(pReal), dimension(3,3,3,3), intent(out) :: dLp_dMp
2018-11-30 12:16:26 +05:30
integer(pInt) :: instance,of,i,k,l,m,n
2018-11-28 00:19:04 +05:30
2015-01-21 20:44:00 +05:30
real(pReal), dimension(plastic_disloUCLA_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))) :: &
2018-11-27 23:31:55 +05:30
gdot_slip_pos,gdot_slip_neg,tau_slip_pos,tau_slip_neg,dgdot_dtauslip_pos,dgdot_dtauslip_neg
!* Shortened notation
of = phasememberAt(ipc,ip,el)
2018-11-30 12:16:26 +05:30
instance = phase_plasticityInstance(phaseAt(ipc,ip,el))
2018-11-29 11:57:35 +05:30
associate(prm => param(instance))
Lp = 0.0_pReal
2018-11-28 11:14:32 +05:30
dLp_dMp = 0.0_pReal
2018-11-30 12:16:26 +05:30
call kinetics(Mp,Temperature,instance,of, &
2018-11-27 23:49:36 +05:30
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos,gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg)
slipSystems: do i = 1_pInt, prm%totalNslip
Lp = Lp + (gdot_slip_pos(i)+gdot_slip_neg(i))*prm%Schmid_slip(1:3,1:3,i)
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
+ dgdot_dtauslip_pos(i) * prm%Schmid_slip(k,l,i) * prm%nonSchmid_pos(m,n,i) &
+ dgdot_dtauslip_neg(i) * prm%Schmid_slip(k,l,i) * prm%nonSchmid_neg(m,n,i)
enddo slipSystems
2018-11-29 02:52:13 +05:30
end associate
Lp = 0.5_pReal * Lp
dLp_dMp = 0.5_pReal * dLp_dMp
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: &
pi
implicit none
2018-11-28 11:14:32 +05:30
real(pReal), dimension(3,3), intent(in):: &
Mp !< 2nd Piola Kirchhoff stress tensor in Mandel notation
real(pReal), intent(in) :: &
temperature !< temperature at integration point
integer(pInt), intent(in) :: &
2018-11-30 14:34:41 +05:30
instance, of
integer(pInt) :: ns,j
real(pReal) :: &
EdgeDipMinDistance,&
AtomicVolume,&
VacancyDiffusion,&
DotRhoMultiplication,&
EdgeDipDistance, &
DotRhoEdgeDipAnnihilation, &
DotRhoEdgeEdgeAnnihilation, &
ClimbVelocity, &
DotRhoEdgeDipClimb, &
2018-11-28 00:19:04 +05:30
DotRhoDipFormation
2018-11-30 14:34:41 +05:30
real(pReal), dimension(plastic_disloUCLA_totalNslip(instance)) :: &
2018-11-27 23:11:33 +05:30
gdot_slip_pos, gdot_slip_neg,&
tau_slip_pos,&
2018-11-27 23:23:01 +05:30
tau_slip_neg, &
dgdot_dtauslip_neg,dgdot_dtauslip_pos
2016-04-25 02:06:35 +05:30
ns = plastic_disloUCLA_totalNslip(instance)
2018-11-30 14:34:41 +05:30
dotState(instance)%whole(:,of) = 0.0_pReal
2016-04-25 02:06:35 +05:30
associate(prm => param(instance), stt => state(instance),mse => microstructure(instance))
!* Dislocation density evolution
2018-11-30 12:16:26 +05:30
call kinetics(Mp,Temperature,instance,of, &
2018-11-27 23:23:01 +05:30
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos,gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg)
dotState(instance)%accshear_slip(:,of) = (gdot_slip_pos+gdot_slip_neg)*0.5_pReal
2018-11-30 12:55:23 +05:30
do j = 1_pInt, prm%totalNslip
2018-11-28 00:19:04 +05:30
!* Multiplication
2018-11-28 00:19:04 +05:30
DotRhoMultiplication = abs(dotState(instance)%accshear_slip(j,of))/&
(prm%burgers(j)* &
mse%mfp(j,of))
!* Dipole formation
EdgeDipMinDistance = &
plastic_disloUCLA_CEdgeDipMinDistance(instance)*prm%burgers(j)
2018-11-27 23:11:33 +05:30
if (dEq0(tau_slip_pos(j))) then
DotRhoDipFormation = 0.0_pReal
else
EdgeDipDistance = &
2018-11-30 12:55:23 +05:30
(3.0_pReal*prm%mu*prm%burgers(j))/&
2018-11-27 23:11:33 +05:30
(16.0_pReal*pi*abs(tau_slip_pos(j)))
if (EdgeDipDistance>mse%mfp(j,of)) EdgeDipDistance=mse%mfp(j,of)
if (EdgeDipDistance<EdgeDipMinDistance) EdgeDipDistance=EdgeDipMinDistance
DotRhoDipFormation = &
((2.0_pReal*EdgeDipDistance)/prm%burgers(j))*&
2018-11-29 11:57:35 +05:30
stt%rhoEdge(j,of)*abs(dotState(instance)%accshear_slip(j,of))*plastic_disloUCLA_dipoleFormationFactor(instance)
endif
!* Spontaneous annihilation of 2 single edge dislocations
DotRhoEdgeEdgeAnnihilation = &
((2.0_pReal*EdgeDipMinDistance)/prm%burgers(j))*&
2018-11-29 11:57:35 +05:30
stt%rhoEdge(j,of)*abs(dotState(instance)%accshear_slip(j,of))
!* Spontaneous annihilation of a single edge dislocation with a dipole constituent
DotRhoEdgeDipAnnihilation = &
((2.0_pReal*EdgeDipMinDistance)/prm%burgers(j))*&
2018-11-29 11:57:35 +05:30
stt%rhoEdgeDip(j,of)*abs(dotState(instance)%accshear_slip(j,of))
!* Dislocation dipole climb
AtomicVolume = &
plastic_disloUCLA_CAtomicVolume(instance)*prm%burgers(j)**(3.0_pReal)
2018-12-03 15:55:29 +05:30
VacancyDiffusion = prm%D0*exp(-plastic_disloUCLA_Qsd(instance)/(kB*Temperature))
2018-11-27 23:11:33 +05:30
if (dEq0(tau_slip_pos(j))) then
DotRhoEdgeDipClimb = 0.0_pReal
else
ClimbVelocity = &
2018-11-30 12:55:23 +05:30
((3.0_pReal*prm%mu*VacancyDiffusion*AtomicVolume)/(2.0_pReal*pi*kB*Temperature))*&
(1/(EdgeDipDistance+EdgeDipMinDistance))
DotRhoEdgeDipClimb = &
2018-11-29 11:57:35 +05:30
(4.0_pReal*ClimbVelocity*stt%rhoEdgeDip(j,of))/(EdgeDipDistance-EdgeDipMinDistance)
endif
!* Edge dislocation density rate of change
dotState(instance)%rhoEdge(j,of) = &
DotRhoMultiplication-DotRhoDipFormation-DotRhoEdgeEdgeAnnihilation
!* Edge dislocation dipole density rate of change
dotState(instance)%rhoEdgeDip(j,of) = &
DotRhoDipFormation-DotRhoEdgeDipAnnihilation-DotRhoEdgeDipClimb
2018-11-28 00:19:04 +05:30
2018-11-30 12:55:23 +05:30
enddo
end associate
2015-01-21 20:44:00 +05:30
end subroutine plastic_disloUCLA_dotState
!--------------------------------------------------------------------------------------------------
!> @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: &
2018-11-28 11:14:32 +05:30
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) :: &
Temperature !< Mandel stress
integer(pInt), intent(in) :: &
instance, &
of
2018-12-03 15:18:37 +05:30
real(pReal), dimension(sum(plastic_disloUCLA_sizePostResult(:,instance))) :: &
postResults
integer(pInt) :: &
2018-12-03 15:18:37 +05:30
o,c,i
real(pReal), dimension(param(instance)%totalNslip) :: &
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos, &
gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg
2018-11-30 14:34:41 +05:30
2018-12-03 15:18:37 +05:30
associate( prm => param(instance), stt => state(instance), mse => microstructure(instance))
2018-11-30 14:34:41 +05:30
postResults = 0.0_pReal
2018-12-03 15:18:37 +05:30
c = 0_pInt
outputsLoop: do o = 1_pInt,size(prm%outputID)
select case(prm%outputID(o))
case (rho_ID)
postResults(c+1_pInt:c+prm%totalNslip) = stt%rhoEdge(1_pInt:prm%totalNslip,of)
case (rhoDip_ID)
postResults(c+1_pInt:c+prm%totalNslip) = stt%rhoEdgeDip(1_pInt:prm%totalNslip,of)
case (shearrate_ID,stressexponent_ID)
call kinetics(Mp,Temperature,instance,of, &
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos,gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg)
if (prm%outputID(o) == shearrate_ID) then
postResults(c+1:c+prm%totalNslip) = (gdot_slip_pos + gdot_slip_neg)*0.5_pReal
elseif(prm%outputID(o) == stressexponent_ID) then
where (dNeq0(gdot_slip_pos+gdot_slip_neg))
postResults(c+1_pInt:c + prm%totalNslip) = (tau_slip_pos+tau_slip_neg) * 0.5_pReal &
/ (gdot_slip_pos+gdot_slip_neg) &
* (dgdot_dtauslip_pos+dgdot_dtauslip_neg)
else where
postResults(c+1_pInt:c + prm%totalNslip) = 0.0_pReal
end where
endif
case (accumulatedshear_ID)
postResults(c+1_pInt:c+prm%totalNslip) = stt%accshear_slip(1_pInt:prm%totalNslip, of)
case (mfp_ID)
postResults(c+1_pInt:c+prm%totalNslip) = mse%mfp(1_pInt:prm%totalNslip, of)
case (resolvedstress_ID)
do i = 1_pInt, prm%totalNslip
postResults(c+i) =math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i))
enddo
case (thresholdstress_ID)
postResults(c+1_pInt:c+prm%totalNslip) = mse%threshold_stress(1_pInt:prm%totalNslip,of)
case (dipoleDistance_ID)
do i = 1_pInt, prm%totalNslip
if (dNeq0(abs(math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i))))) then
postResults(c+i) = (3.0_pReal*prm%mu*prm%burgers(i)) &
/ (16.0_pReal*pi*abs(math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i))))
else
postResults(c+i) = huge(1.0_pReal)
endif
postResults(c+i)=min(postResults(c+i),mse%mfp(i,of))
enddo
end select
c = c + prm%totalNslip
enddo outputsLoop
end associate
end function plastic_disloUCLA_postResults
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2018-11-30 12:16:26 +05:30
subroutine kinetics(Mp,Temperature,instance,of, &
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos,gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg)
use prec, only: &
tol_math_check, &
dEq, dNeq0
use math, only: &
2018-11-28 11:14:32 +05:30
pi, &
math_mul33xx33
implicit none
2018-11-28 11:14:32 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< 2nd Piola Kirchhoff stress tensor in Mandel notation
real(pReal), intent(in) :: &
temperature !< temperature at integration point
integer(pInt), intent(in) :: &
2018-11-30 12:16:26 +05:30
instance,of
integer(pInt) :: &
2018-11-30 12:55:23 +05:30
i,j
real(pReal) :: StressRatio_p,StressRatio_pminus1,&
BoltzmannRatio,DotGamma0,stressRatio,&
dvel_slip, vel_slip
real(pReal), intent(out), dimension(plastic_disloUCLA_totalNslip(instance)) :: &
gdot_slip_pos,dgdot_dtauslip_pos,tau_slip_pos,gdot_slip_neg,dgdot_dtauslip_neg,tau_slip_neg
associate(prm => param(instance), stt => state(instance),mse => microstructure(instance))
2018-11-30 12:16:26 +05:30
gdot_slip_pos = 0.0_pReal
gdot_slip_neg = 0.0_pReal
dgdot_dtauslip_pos = 0.0_pReal
dgdot_dtauslip_neg = 0.0_pReal
2018-11-30 12:16:26 +05:30
do j = 1_pInt, prm%totalNslip
2014-11-07 16:53:34 +05:30
!* Boltzmann ratio
BoltzmannRatio = prm%H0kp(j)/(kB*Temperature)
2014-11-07 16:53:34 +05:30
!* Initial shear rates
2018-11-29 15:01:02 +05:30
DotGamma0 = stt%rhoEdge(j,of)*prm%burgers(j)*prm%v0(j)
!* Resolved shear stress on slip system
tau_slip_pos(j) = math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,j))
tau_slip_neg(j) = math_mul33xx33(Mp,prm%nonSchmid_neg(1:3,1:3,j))
2014-11-07 16:53:34 +05:30
significantPositiveTau: if((abs(tau_slip_pos(j))-mse%threshold_stress(j, of)) > tol_math_check) then
2015-01-21 20:44:00 +05:30
!* Stress ratio
stressRatio = ((abs(tau_slip_pos(j))-mse%threshold_stress(j, of))/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+&
2018-11-29 13:14:31 +05:30
prm%tau_Peierls(j)))
stressRatio_p = stressRatio** prm%p(j)
stressRatio_pminus1 = stressRatio**(prm%p(j)-1.0_pReal)
2015-01-21 20:44:00 +05:30
!* Shear rates due to slip
vel_slip = 2.0_pReal*prm%burgers(j) &
* prm%kink_height(j) * prm%omega(j) &
* ( mse%mfp(j,of) - prm%kink_width(j) ) &
2015-01-16 04:10:17 +05:30
* (tau_slip_pos(j) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) ) &
2015-01-16 04:10:17 +05:30
/ ( &
2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_pos(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
)
gdot_slip_pos(j) = DotGamma0 &
2015-01-15 16:26:15 +05:30
* vel_slip &
* sign(1.0_pReal,tau_slip_pos(j))
!* Derivatives of shear rates
2015-01-21 20:44:00 +05:30
dvel_slip = &
2.0_pReal*prm%burgers(j) &
* prm%kink_height(j) * prm%omega(j) &
* ( mse%mfp(j,of) - prm%kink_width(j) ) &
2015-01-16 04:10:17 +05:30
* ( &
(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
+ tau_slip_pos(j) &
* (abs(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)))&
*BoltzmannRatio*prm%p(j)&
*prm%q(j)/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+prm%tau_Peierls(j))*&
StressRatio_pminus1*(1-StressRatio_p)**(prm%q(j)-1.0_pReal) ) &
2015-01-16 04:10:17 +05:30
) &
* (2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_pos(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
) &
- (tau_slip_pos(j) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) ) &
* (2.0_pReal*(prm%burgers(j)**2.0_pReal) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* (abs(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)))&
*BoltzmannRatio*prm%p(j)&
*prm%q(j)/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+prm%tau_Peierls(j))*&
StressRatio_pminus1*(1-StressRatio_p)**(prm%q(j)-1.0_pReal) )&
2015-01-16 04:10:17 +05:30
) &
) &
/ ( &
( &
2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_pos(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
)**2.0_pReal &
)
2015-01-15 16:26:15 +05:30
dgdot_dtauslip_pos(j) = DotGamma0 * dvel_slip
endif significantPositiveTau
2018-11-30 12:16:26 +05:30
significantNegativeTau: if((abs(tau_slip_neg(j))-mse%threshold_stress(j, of)) > tol_math_check) then
!* Stress ratios
stressRatio = ((abs(tau_slip_neg(j))-mse%threshold_stress(j, of))/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+&
2018-11-29 13:14:31 +05:30
prm%tau_Peierls(j)))
stressRatio_p = stressRatio** prm%p(j)
stressRatio_pminus1 = stressRatio**(prm%p(j)-1.0_pReal)
!* Shear rates due to slip
vel_slip = 2.0_pReal*prm%burgers(j) &
* prm%kink_height(j) * prm%omega(j) &
* ( mse%mfp(j,of) - prm%kink_width(j) ) &
2015-01-16 04:10:17 +05:30
* (tau_slip_neg(j) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) ) &
2015-01-16 04:10:17 +05:30
/ ( &
2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_neg(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
)
2015-01-15 16:26:15 +05:30
gdot_slip_neg(j) = DotGamma0 &
2015-01-15 16:26:15 +05:30
* vel_slip &
* sign(1.0_pReal,tau_slip_neg(j))
!* Derivatives of shear rates
dvel_slip = &
2.0_pReal*prm%burgers(j) &
* prm%kink_height(j) * prm%omega(j) &
* ( mse%mfp(j,of) - prm%kink_width(j) ) &
2015-01-16 04:10:17 +05:30
* ( &
(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
+ tau_slip_neg(j) &
* (abs(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)))&
*BoltzmannRatio*prm%p(j)&
*prm%q(j)/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+prm%tau_Peierls(j))*&
StressRatio_pminus1*(1-StressRatio_p)**(prm%q(j)-1.0_pReal) ) &
2015-01-16 04:10:17 +05:30
) &
* (2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_neg(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
) &
- (tau_slip_neg(j) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) ) &
* (2.0_pReal*(prm%burgers(j)**2.0_pReal) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* (abs(exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)))&
*BoltzmannRatio*prm%p(j)&
*prm%q(j)/&
2018-11-29 15:01:02 +05:30
(prm%solidSolutionStrength+prm%tau_Peierls(j))*&
StressRatio_pminus1*(1-StressRatio_p)**(prm%q(j)-1.0_pReal) )&
2015-01-16 04:10:17 +05:30
) &
) &
/ ( &
( &
2.0_pReal*(prm%burgers(j)**2.0_pReal)*tau_slip_neg(j) &
2018-11-29 15:07:06 +05:30
+ prm%omega(j) * prm%B(j) &
*(( mse%mfp(j,of) - prm%kink_width(j) )**2.0_pReal) &
* exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q(j)) &
2015-01-16 04:10:17 +05:30
)**2.0_pReal &
)
2015-01-15 16:26:15 +05:30
dgdot_dtauslip_neg(j) = DotGamma0 * dvel_slip
2018-11-30 12:16:26 +05:30
endif significantNegativeTau
enddo
end associate
end subroutine kinetics
2015-01-21 20:44:00 +05:30
end module plastic_disloUCLA