DAMASK_EICMD/src/plastic_kinematichardening.f90

622 lines
26 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Philip Eisenlohr, Michigan State University
!> @author Zhuowen Zhao, Michigan State University
2018-11-25 15:44:09 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Phenomenological crystal plasticity using a power law formulation for the shear rates
!! and a Voce-type kinematic hardening rule
!--------------------------------------------------------------------------------------------------
module plastic_kinehardening
use prec
use debug
use math
use IO
use material
use config
use lattice
use discretization
use results
implicit none
private
2019-03-28 02:36:40 +05:30
integer, dimension(:,:), allocatable, target, public :: &
plastic_kinehardening_sizePostResult !< size of each post result output
2019-03-28 02:36:40 +05:30
character(len=64), dimension(:,:), allocatable, target, public :: &
plastic_kinehardening_output !< name of each post result output
enum, bind(c)
2018-11-25 15:44:09 +05:30
enumerator :: &
undefined_ID, &
crss_ID, & !< critical resolved stress
crss_back_ID, & !< critical resolved back stress
sense_ID, & !< sense of acting shear stress (-1 or +1)
chi0_ID, & !< backstress at last switch of stress sense (positive?)
gamma0_ID, & !< accumulated shear at last switch of stress sense (at current switch?)
2018-11-25 15:44:09 +05:30
accshear_ID, &
shearrate_ID, &
resolvedstress_ID
end enum
2018-11-25 15:44:09 +05:30
type :: tParameters
real(pReal) :: &
gdot0, & !< reference shear strain rate for slip
n, & !< stress exponent for slip
aTolResistance, &
aTolShear
2019-03-28 02:36:40 +05:30
real(pReal), allocatable, dimension(:) :: &
crss0, & !< initial critical shear stress for slip
theta0, & !< initial hardening rate of forward stress for each slip
theta1, & !< asymptotic hardening rate of forward stress for each slip
theta0_b, & !< initial hardening rate of back stress for each slip
theta1_b, & !< asymptotic hardening rate of back stress for each slip
tau1, &
tau1_b, &
nonSchmidCoeff
2019-03-28 02:36:40 +05:30
real(pReal), allocatable, dimension(:,:) :: &
interaction_slipslip !< slip resistance from slip activity
2019-03-28 02:36:40 +05:30
real(pReal), allocatable, dimension(:,:,:) :: &
Schmid, &
nonSchmid_pos, &
nonSchmid_neg
2019-03-28 02:36:40 +05:30
integer :: &
2018-12-30 20:09:48 +05:30
totalNslip, & !< total number of active slip system
2019-03-28 02:36:40 +05:30
of_debug = 0
integer, 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 tParameters
type :: tKinehardeningState
real(pReal), pointer, dimension(:,:) :: & !< vectors along NipcMyInstance
crss, & !< critical resolved stress
crss_back, & !< critical resolved back stress
sense, & !< sense of acting shear stress (-1 or +1)
chi0, & !< backstress at last switch of stress sense
gamma0, & !< accumulated shear at last switch of stress sense
accshear !< accumulated (absolute) shear
end type tKinehardeningState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:) :: param
type(tKinehardeningState), allocatable, dimension(:) :: &
dotState, &
deltaState, &
2018-12-22 18:41:34 +05:30
state
public :: &
plastic_kinehardening_init, &
plastic_kinehardening_LpAndItsTangent, &
plastic_kinehardening_dotState, &
plastic_kinehardening_deltaState, &
plastic_kinehardening_postResults, &
plastic_kinehardening_results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2018-12-13 14:21:43 +05:30
subroutine plastic_kinehardening_init
2019-03-28 02:36:40 +05:30
integer :: &
2018-12-22 03:37:31 +05:30
Ninstance, &
p, i, o, &
2019-01-07 12:34:02 +05:30
NipcMyPhase, &
sizeState, sizeDeltaState, sizeDotState, &
startIndex, endIndex
2019-03-28 02:36:40 +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)::]
integer(kind(undefined_ID)) :: &
outputID
character(len=pStringLen) :: &
extmsg = ''
character(len=65536), dimension(:), allocatable :: &
outputs
write(6,'(/,a)') ' <<<+- plastic_'//PLASTICITY_KINEHARDENING_label//' init -+>>>'
2019-03-09 05:37:26 +05:30
Ninstance = count(phase_plasticity == PLASTICITY_KINEHARDENING_ID)
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-03-28 02:36:40 +05:30
allocate(plastic_kinehardening_sizePostResult(maxval(phase_Noutput),Ninstance),source=0)
2018-12-22 03:37:31 +05:30
allocate(plastic_kinehardening_output(maxval(phase_Noutput),Ninstance))
plastic_kinehardening_output = ''
allocate(param(Ninstance))
2018-12-22 03:37:31 +05:30
allocate(state(Ninstance))
allocate(dotState(Ninstance))
allocate(deltaState(Ninstance))
2019-03-28 02:36:40 +05:30
do p = 1, size(phase_plasticityInstance)
if (phase_plasticity(p) /= PLASTICITY_KINEHARDENING_ID) cycle
2018-12-22 19:52:41 +05:30
associate(prm => param(phase_plasticityInstance(p)), &
dot => dotState(phase_plasticityInstance(p)), &
2019-01-07 11:37:55 +05:30
dlt => deltaState(phase_plasticityInstance(p)), &
stt => state(phase_plasticityInstance(p)),&
config => config_phase(p))
2018-12-30 20:09:48 +05:30
#ifdef DEBUG
if (p==material_phaseAt(debug_g,debug_e)) then
2019-06-14 12:47:05 +05:30
prm%of_debug = material_phasememberAt(debug_g,debug_i,debug_e)
2018-12-30 20:09:48 +05:30
endif
#endif
!--------------------------------------------------------------------------------------------------
! optional parameters that need to be defined
prm%aTolResistance = config%getFloat('atol_resistance',defaultVal=1.0_pReal)
prm%aTolShear = config%getFloat('atol_shear', defaultVal=1.0e-6_pReal)
! sanity checks
2019-01-07 11:37:55 +05:30
if (prm%aTolResistance <= 0.0_pReal) extmsg = trim(extmsg)//' aTolresistance'
if (prm%aTolShear <= 0.0_pReal) extmsg = trim(extmsg)//' aTolShear'
!--------------------------------------------------------------------------------------------------
! slip related parameters
prm%Nslip = config%getInts('nslip',defaultVal=emptyIntArray)
prm%totalNslip = sum(prm%Nslip)
2019-03-28 02:36:40 +05:30
slipActive: if (prm%totalNslip > 0) then
prm%Schmid = lattice_SchmidMatrix_slip(prm%Nslip,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',&
defaultVal = emptyRealArray)
2019-03-28 02:36:40 +05:30
prm%nonSchmid_pos = lattice_nonSchmidMatrix(prm%Nslip,prm%nonSchmidCoeff,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(prm%Nslip,prm%nonSchmidCoeff,-1)
else
prm%nonSchmid_pos = prm%Schmid
prm%nonSchmid_neg = prm%Schmid
endif
prm%interaction_SlipSlip = lattice_interaction_SlipBySlip(prm%Nslip, &
config%getFloats('interaction_slipslip'), &
config%getString('lattice_structure'))
prm%crss0 = config%getFloats('crss0', requiredSize=size(prm%Nslip))
prm%tau1 = config%getFloats('tau1', requiredSize=size(prm%Nslip))
prm%tau1_b = config%getFloats('tau1_b', requiredSize=size(prm%Nslip))
prm%theta0 = config%getFloats('theta0', requiredSize=size(prm%Nslip))
prm%theta1 = config%getFloats('theta1', requiredSize=size(prm%Nslip))
prm%theta0_b = config%getFloats('theta0_b', requiredSize=size(prm%Nslip))
prm%theta1_b = config%getFloats('theta1_b', requiredSize=size(prm%Nslip))
2019-01-07 11:54:02 +05:30
2019-01-07 11:37:55 +05:30
prm%gdot0 = config%getFloat('gdot0')
prm%n = config%getFloat('n_slip')
2018-12-13 12:01:56 +05:30
! expand: family => system
2019-01-07 11:37:55 +05:30
prm%crss0 = math_expand(prm%crss0, prm%Nslip)
prm%tau1 = math_expand(prm%tau1, prm%Nslip)
prm%tau1_b = math_expand(prm%tau1_b, prm%Nslip)
prm%theta0 = math_expand(prm%theta0, prm%Nslip)
prm%theta1 = math_expand(prm%theta1, prm%Nslip)
2018-12-13 12:01:56 +05:30
prm%theta0_b = math_expand(prm%theta0_b,prm%Nslip)
prm%theta1_b = math_expand(prm%theta1_b,prm%Nslip)
2019-01-07 11:37:55 +05:30
!--------------------------------------------------------------------------------------------------
! sanity checks
if ( prm%gdot0 <= 0.0_pReal) extmsg = trim(extmsg)//' gdot0'
if ( prm%n <= 0.0_pReal) extmsg = trim(extmsg)//' n_slip'
if (any(prm%crss0 <= 0.0_pReal)) extmsg = trim(extmsg)//' crss0'
2019-01-07 12:34:02 +05:30
if (any(prm%tau1 <= 0.0_pReal)) extmsg = trim(extmsg)//' tau1'
if (any(prm%tau1_b <= 0.0_pReal)) extmsg = trim(extmsg)//' tau1_b'
!ToDo: Any sensible checks for theta?
2018-11-26 11:40:43 +05:30
endif slipActive
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') &
2019-03-28 02:36:40 +05:30
call IO_error(211,ext_msg=trim(extmsg)//'('//PLASTICITY_KINEHARDENING_label//')')
!--------------------------------------------------------------------------------------------------
! output pararameters
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
allocate(prm%outputID(0))
2019-03-28 02:36:40 +05:30
do i=1, size(outputs)
outputID = undefined_ID
select case(outputs(i))
2019-01-07 11:54:02 +05:30
2019-01-07 12:34:02 +05:30
case ('resistance')
2019-03-28 02:36:40 +05:30
outputID = merge(crss_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('accumulatedshear')
2019-03-28 02:36:40 +05:30
outputID = merge(accshear_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('shearrate')
2019-03-28 02:36:40 +05:30
outputID = merge(shearrate_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('resolvedstress')
2019-03-28 02:36:40 +05:30
outputID = merge(resolvedstress_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('backstress')
2019-03-28 02:36:40 +05:30
outputID = merge(crss_back_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('sense')
2019-03-28 02:36:40 +05:30
outputID = merge(sense_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('chi0')
2019-03-28 02:36:40 +05:30
outputID = merge(chi0_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
case ('gamma0')
2019-03-28 02:36:40 +05:30
outputID = merge(gamma0_ID,undefined_ID,prm%totalNslip>0)
2019-01-07 12:34:02 +05:30
end select
if (outputID /= undefined_ID) then
plastic_kinehardening_output(i,phase_plasticityInstance(p)) = outputs(i)
plastic_kinehardening_sizePostResult(i,phase_plasticityInstance(p)) = prm%totalNslip
prm%outputID = [prm%outputID , outputID]
endif
enddo
2018-12-22 19:52:41 +05:30
!--------------------------------------------------------------------------------------------------
! allocate state arrays
NipcMyPhase = count(material_phaseAt == p) * discretization_nIP
2019-01-07 12:34:02 +05:30
sizeDotState = size(['crss ','crss_back', 'accshear ']) * prm%totalNslip
sizeDeltaState = size(['sense ', 'chi0 ', 'gamma0' ]) * prm%totalNslip
2018-12-22 03:37:31 +05:30
sizeState = sizeDotState + sizeDeltaState
call material_allocatePlasticState(p,NipcMyPhase,sizeState,sizeDotState,sizeDeltaState, &
2019-03-28 02:36:40 +05:30
prm%totalNslip,0,0)
2018-12-22 03:37:31 +05:30
plasticState(p)%sizePostResults = sum(plastic_kinehardening_sizePostResult(:,phase_plasticityInstance(p)))
!--------------------------------------------------------------------------------------------------
! locally defined state aliases and initialization of state0 and aTolState
2019-03-28 02:36:40 +05:30
startIndex = 1
endIndex = prm%totalNslip
stt%crss => plasticState(p)%state (startIndex:endIndex,:)
stt%crss = spread(prm%crss0, 2, NipcMyPhase)
dot%crss => plasticState(p)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolResistance
2019-03-28 02:36:40 +05:30
startIndex = endIndex + 1
2019-01-07 11:54:02 +05:30
endIndex = endIndex + prm%totalNslip
stt%crss_back => plasticState(p)%state (startIndex:endIndex,:)
dot%crss_back => plasticState(p)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolResistance
2019-03-28 02:36:40 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%totalNslip
stt%accshear => plasticState(p)%state (startIndex:endIndex,:)
dot%accshear => plasticState(p)%dotState(startIndex:endIndex,:)
plasticState(p)%aTolState(startIndex:endIndex) = prm%aTolShear
! global alias
plasticState(p)%slipRate => plasticState(p)%dotState(startIndex:endIndex,:)
plasticState(p)%accumulatedSlip => plasticState(p)%state(startIndex:endIndex,:)
o = plasticState(p)%offsetDeltaState
2019-03-28 02:36:40 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%totalNslip
stt%sense => plasticState(p)%state (startIndex :endIndex ,:)
dlt%sense => plasticState(p)%deltaState(startIndex-o:endIndex-o,:)
2019-01-07 11:54:02 +05:30
2019-03-28 02:36:40 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%totalNslip
stt%chi0 => plasticState(p)%state (startIndex :endIndex ,:)
dlt%chi0 => plasticState(p)%deltaState(startIndex-o:endIndex-o,:)
2019-03-28 02:36:40 +05:30
startIndex = endIndex + 1
endIndex = endIndex + prm%totalNslip
stt%gamma0 => plasticState(p)%state (startIndex :endIndex ,:)
dlt%gamma0 => plasticState(p)%deltaState(startIndex-o:endIndex-o,:)
2019-01-07 11:54:02 +05:30
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
2018-12-22 03:37:31 +05:30
end associate
2019-01-07 12:34:02 +05:30
enddo
end subroutine plastic_kinehardening_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
!--------------------------------------------------------------------------------------------------
2018-12-25 18:50:01 +05:30
pure subroutine plastic_kinehardening_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
real(pReal), dimension(3,3), intent(out) :: &
Lp !< plastic velocity gradient
2018-11-23 10:07:31 +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-13 14:04:40 +05:30
Mp !< Mandel stress
2019-03-28 02:36:40 +05:30
integer, intent(in) :: &
2018-12-13 14:04:40 +05:30
instance, &
of
2019-03-28 02:36:40 +05:30
integer :: &
i,k,l,m,n
2018-12-22 19:52:41 +05:30
real(pReal), dimension(param(instance)%totalNslip) :: &
gdot_pos,gdot_neg, &
dgdot_dtau_pos,dgdot_dtau_neg
Lp = 0.0_pReal
2018-11-23 10:07:31 +05:30
dLp_dMp = 0.0_pReal
2019-01-07 12:34:02 +05:30
associate(prm => param(instance))
2018-12-22 19:52:41 +05:30
call kinetics(Mp,instance,of,gdot_pos,gdot_neg,dgdot_dtau_pos,dgdot_dtau_neg)
2019-03-28 02:36:40 +05:30
do i = 1, prm%totalNslip
Lp = Lp + (gdot_pos(i)+gdot_neg(i))*prm%Schmid(1:3,1:3,i)
2019-03-28 02:36:40 +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)
2018-12-13 12:01:56 +05:30
enddo
2019-01-07 11:54:02 +05:30
2018-12-25 18:50:01 +05:30
end associate
end subroutine plastic_kinehardening_LpAndItsTangent
2018-12-22 03:11:39 +05:30
2019-01-07 12:06:11 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates the rate of change of microstructure
!--------------------------------------------------------------------------------------------------
subroutine plastic_kinehardening_dotState(Mp,instance,of)
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2019-03-28 02:36:40 +05:30
integer, intent(in) :: &
2019-01-07 12:06:11 +05:30
instance, &
of
real(pReal) :: &
sumGamma
2019-01-07 12:34:02 +05:30
real(pReal), dimension(param(instance)%totalNslip) :: &
gdot_pos,gdot_neg
2019-01-07 12:06:11 +05:30
associate(prm => param(instance), stt => state(instance), dot => dotState(instance))
call kinetics(Mp,instance,of,gdot_pos,gdot_neg)
dot%accshear(:,of) = abs(gdot_pos+gdot_neg)
sumGamma = sum(stt%accshear(:,of))
2019-03-26 11:17:17 +05:30
dot%crss(:,of) = matmul(prm%interaction_SlipSlip,dot%accshear(:,of)) &
* ( prm%theta1 &
+ (prm%theta0 - prm%theta1 + prm%theta0*prm%theta1*sumGamma/prm%tau1) &
* exp(-sumGamma*prm%theta0/prm%tau1) &
)
2019-01-07 12:06:11 +05:30
dot%crss_back(:,of) = stt%sense(:,of)*dot%accshear(:,of) * &
( prm%theta1_b + &
(prm%theta0_b - prm%theta1_b &
+ prm%theta0_b*prm%theta1_b/(prm%tau1_b+stt%chi0(:,of))*(stt%accshear(:,of)-stt%gamma0(:,of))&
) *exp(-(stt%accshear(:,of)-stt%gamma0(:,of)) *prm%theta0_b/(prm%tau1_b+stt%chi0(:,of))) &
)
end associate
end subroutine plastic_kinehardening_dotState
!--------------------------------------------------------------------------------------------------
!> @brief calculates (instantaneous) incremental change of microstructure
!--------------------------------------------------------------------------------------------------
2018-12-13 14:04:40 +05:30
subroutine plastic_kinehardening_deltaState(Mp,instance,of)
2018-12-30 20:09:48 +05:30
2018-12-13 14:04:40 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2019-03-28 02:36:40 +05:30
integer, intent(in) :: &
2018-12-13 14:04:40 +05:30
instance, &
of
2018-12-22 19:52:41 +05:30
real(pReal), dimension(param(instance)%totalNslip) :: &
gdot_pos,gdot_neg, &
sense
2019-01-07 11:37:55 +05:30
associate(prm => param(instance), stt => state(instance), dlt => deltaState(instance))
2019-01-07 11:54:02 +05:30
2018-12-22 19:52:41 +05:30
call kinetics(Mp,instance,of,gdot_pos,gdot_neg)
sense = merge(state(instance)%sense(:,of), & ! keep existing...
2019-01-07 11:54:02 +05:30
sign(1.0_pReal,gdot_pos+gdot_neg), & ! ...or have a defined
dEq0(gdot_pos+gdot_neg,1e-10_pReal)) ! current sense of shear direction
2017-11-07 04:41:02 +05:30
#ifdef DEBUG
2019-03-28 02:36:40 +05:30
if (iand(debug_level(debug_constitutive), debug_levelExtensive) /= 0 &
2018-12-30 20:09:48 +05:30
.and. (of == prm%of_debug &
2019-03-28 02:36:40 +05:30
.or. .not. iand(debug_level(debug_constitutive),debug_levelSelective) /= 0)) then
2018-12-30 20:09:48 +05:30
write(6,'(a)') '======= kinehardening delta state ======='
write(6,*) sense,state(instance)%sense(:,of)
endif
2017-11-07 04:41:02 +05:30
#endif
2018-12-22 19:52:41 +05:30
2018-12-13 12:01:56 +05:30
!--------------------------------------------------------------------------------------------------
! switch in sense of shear?
2018-12-22 19:17:02 +05:30
where(dNeq(sense,stt%sense(:,of),0.1_pReal))
2019-01-07 11:37:55 +05:30
dlt%sense (:,of) = sense - stt%sense(:,of) ! switch sense
dlt%chi0 (:,of) = abs(stt%crss_back(:,of)) - stt%chi0(:,of) ! remember current backstress magnitude
dlt%gamma0(:,of) = stt%accshear(:,of) - stt%gamma0(:,of) ! remember current accumulated shear
2018-12-13 12:01:56 +05:30
else where
2019-01-07 11:37:55 +05:30
dlt%sense (:,of) = 0.0_pReal
dlt%chi0 (:,of) = 0.0_pReal
dlt%gamma0(:,of) = 0.0_pReal
2018-12-13 12:01:56 +05:30
end where
2019-01-07 11:54:02 +05:30
2018-12-22 19:17:02 +05:30
end associate
end subroutine plastic_kinehardening_deltaState
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2018-12-13 14:04:40 +05:30
function plastic_kinehardening_postResults(Mp,instance,of) result(postResults)
real(pReal), dimension(3,3), intent(in) :: &
2018-12-13 14:04:40 +05:30
Mp !< Mandel stress
2019-03-28 02:36:40 +05:30
integer, intent(in) :: &
2018-12-13 14:04:40 +05:30
instance, &
of
2018-12-13 14:04:40 +05:30
real(pReal), dimension(sum(plastic_kinehardening_sizePostResult(:,instance))) :: &
postResults
2019-03-28 02:36:40 +05:30
integer :: &
o,c,i
2018-12-22 19:52:41 +05:30
real(pReal), dimension(param(instance)%totalNslip) :: &
2018-12-22 18:41:34 +05:30
gdot_pos,gdot_neg
2019-03-28 02:36:40 +05:30
c = 0
2018-12-22 18:41:34 +05:30
associate(prm => param(instance), stt => state(instance))
2019-03-28 02:36:40 +05:30
outputsLoop: do o = 1,size(prm%outputID)
2018-12-13 11:29:56 +05:30
select case(prm%outputID(o))
case (crss_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%crss(:,of)
case(crss_back_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%crss_back(:,of)
case (sense_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%sense(:,of)
case (chi0_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%chi0(:,of)
case (gamma0_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%gamma0(:,of)
case (accshear_ID)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = stt%accshear(:,of)
case (shearrate_ID)
2019-01-07 12:34:02 +05:30
call kinetics(Mp,instance,of,gdot_pos,gdot_neg)
2019-03-28 02:36:40 +05:30
postResults(c+1:c+prm%totalNslip) = gdot_pos+gdot_neg
case (resolvedstress_ID)
2019-03-28 02:36:40 +05:30
do i = 1, prm%totalNslip
postResults(c+i) = math_mul33xx33(Mp,prm%Schmid(1:3,1:3,i))
2018-12-13 12:01:56 +05:30
enddo
end select
2019-01-07 11:54:02 +05:30
2019-01-07 11:37:55 +05:30
c = c + prm%totalNslip
2019-01-07 11:54:02 +05:30
enddo outputsLoop
2018-12-13 11:29:56 +05:30
end associate
end function plastic_kinehardening_postResults
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine plastic_kinehardening_results(instance,group)
2019-04-04 13:34:44 +05:30
#if defined(PETSc) || defined(DAMASK_HDF5)
integer, intent(in) :: instance
character(len=*) :: group
integer :: o
associate(prm => param(instance), stt => state(instance))
2019-03-28 02:36:40 +05:30
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
case (crss_ID)
call results_writeDataset(group,stt%crss,'xi_sl', &
'resistance against plastic slip','Pa')
case(crss_back_ID)
call results_writeDataset(group,stt%crss_back,'tau_back', &
'back stress against plastic slip','Pa')
case (sense_ID)
call results_writeDataset(group,stt%sense,'sense_of_shear','tbd','1')
case (chi0_ID)
call results_writeDataset(group,stt%chi0,'chi0','tbd','Pa')
case (gamma0_ID)
call results_writeDataset(group,stt%gamma0,'gamma0','tbd','1')
case (accshear_ID)
call results_writeDataset(group,stt%accshear,'gamma_sl', &
'plastic shear','1')
end select
enddo outputsLoop
end associate
#else
integer, intent(in) :: instance
character(len=*) :: group
#endif
end subroutine plastic_kinehardening_results
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on slip systems and derivatives with respect to resolved stress
!> @details: Shear rates are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end
2018-11-26 11:40:43 +05:30
!--------------------------------------------------------------------------------------------------
pure subroutine kinetics(Mp,instance,of, &
gdot_pos,gdot_neg,dgdot_dtau_pos,dgdot_dtau_neg)
2018-11-26 11:40:43 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
2019-03-28 02:36:40 +05:30
integer, intent(in) :: &
instance, &
2018-11-26 11:40:43 +05:30
of
real(pReal), intent(out), dimension(param(instance)%totalNslip) :: &
2018-11-26 11:40:43 +05:30
gdot_pos, &
gdot_neg
real(pReal), intent(out), optional, dimension(param(instance)%totalNslip) :: &
2018-11-26 11:40:43 +05:30
dgdot_dtau_pos, &
dgdot_dtau_neg
2018-12-22 19:52:41 +05:30
real(pReal), dimension(param(instance)%totalNslip) :: &
2018-11-26 11:40:43 +05:30
tau_pos, &
tau_neg
2019-03-28 02:36:40 +05:30
integer :: i
logical :: nonSchmidActive
2018-11-26 11:40:43 +05:30
associate(prm => param(instance), stt => state(instance))
2019-01-07 11:54:02 +05:30
2019-03-28 02:36:40 +05:30
nonSchmidActive = size(prm%nonSchmidCoeff) > 0
2018-11-26 11:40:43 +05:30
2019-03-28 02:36:40 +05:30
do i = 1, prm%totalNslip
2018-12-22 19:52:41 +05:30
tau_pos(i) = math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i)) - stt%crss_back(i,of)
tau_neg(i) = merge(math_mul33xx33(Mp,prm%nonSchmid_neg(1:3,1:3,i)) - stt%crss_back(i,of), &
2018-12-25 18:50:01 +05:30
0.0_pReal, nonSchmidActive)
2018-11-26 11:40:43 +05:30
enddo
where(dNeq0(tau_pos))
2018-11-27 03:06:32 +05:30
gdot_pos = prm%gdot0 * merge(0.5_pReal,1.0_pReal, nonSchmidActive) & ! 1/2 if non-Schmid active
* sign(abs(tau_pos/stt%crss(:,of))**prm%n, tau_pos)
2018-11-26 11:40:43 +05:30
else where
gdot_pos = 0.0_pReal
end where
where(dNeq0(tau_neg))
2018-11-27 03:06:32 +05:30
gdot_neg = prm%gdot0 * 0.5_pReal & ! only used if non-Schmid active, always 1/2
* sign(abs(tau_neg/stt%crss(:,of))**prm%n, tau_neg)
2018-11-26 11:40:43 +05:30
else where
gdot_neg = 0.0_pReal
end where
if (present(dgdot_dtau_pos)) then
where(dNeq0(gdot_pos))
dgdot_dtau_pos = gdot_pos*prm%n/tau_pos
2018-11-26 11:40:43 +05:30
else where
dgdot_dtau_pos = 0.0_pReal
end where
endif
if (present(dgdot_dtau_neg)) then
where(dNeq0(gdot_neg))
dgdot_dtau_neg = gdot_neg*prm%n/tau_neg
2018-11-26 11:40:43 +05:30
else where
dgdot_dtau_neg = 0.0_pReal
end where
endif
2018-12-22 03:37:31 +05:30
end associate
2018-11-26 11:40:43 +05:30
end subroutine kinetics
end module plastic_kinehardening