DAMASK_EICMD/src/phase_mechanical_plastic_ph...

544 lines
24 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 Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief phenomenological crystal plasticity formulation using a powerlaw fitting
!--------------------------------------------------------------------------------------------------
submodule(phase:plastic) phenopowerlaw
type :: tParameters
real(pREAL) :: &
dot_gamma_0_sl = 1.0_pREAL, & !< reference shear strain rate for slip
dot_gamma_0_tw = 1.0_pREAL, & !< reference shear strain rate for twin
n_sl = 1.0_pREAL, & !< stress exponent for slip
n_tw = 1.0_pREAL, & !< stress exponent for twin
f_sat_sl_tw = 1.0_pREAL, & !< push-up factor for slip saturation due to twinning
c_1 = 1.0_pREAL, &
c_2 = 1.0_pREAL, &
c_3 = 1.0_pREAL, &
c_4 = 1.0_pREAL, &
h_0_sl_sl = 1.0_pREAL, & !< reference hardening slip - slip
h_0_tw_sl = 1.0_pREAL, & !< reference hardening twin - slip
h_0_tw_tw = 1.0_pREAL, & !< reference hardening twin - twin
a_sl = 1.0_pREAL
real(pREAL), allocatable, dimension(:) :: &
2020-09-22 19:34:14 +05:30
xi_inf_sl, & !< maximum critical shear stress for slip
h_int, & !< per family hardening activity (optional)
gamma_char !< characteristic shear for twins
real(pREAL), allocatable, dimension(:,:) :: &
2020-09-22 19:34:14 +05:30
h_sl_sl, & !< slip resistance from slip activity
h_sl_tw, & !< slip resistance from twin activity
h_tw_sl, & !< twin resistance from slip activity
h_tw_tw !< twin resistance from twin activity
real(pREAL), allocatable, dimension(:,:,:) :: &
2020-03-16 01:46:28 +05:30
P_sl, &
P_tw, &
2021-07-21 19:59:57 +05:30
P_nS_pos, &
P_nS_neg
integer :: &
2020-03-16 01:46:28 +05:30
sum_N_sl, & !< total number of active slip system
sum_N_tw !< total number of active twin systems
2020-03-16 16:23:36 +05:30
logical :: &
nonSchmidActive = .false.
2023-06-04 10:47:38 +05:30
character(len=pSTRLEN), allocatable, dimension(:) :: &
output
character(len=:), allocatable, dimension(:) :: &
systems_sl, &
systems_tw
end type tParameters
type :: tIndexDotState
integer, dimension(2) :: &
xi_sl, &
xi_tw, &
gamma_sl, &
gamma_tw
end type tIndexDotState
type :: tPhenopowerlawState
real(pREAL), pointer, dimension(:,:) :: &
2021-07-21 18:33:28 +05:30
xi_sl, &
xi_tw, &
gamma_sl, &
gamma_tw
end type tPhenopowerlawState
!--------------------------------------------------------------------------------------------------
! containers for parameters, dot state index, and state
type(tParameters), allocatable, dimension(:) :: param
type(tIndexDotState), allocatable, dimension(:) :: indexDotState
type(tPhenopowerlawState), allocatable, dimension(:) :: state
contains
!--------------------------------------------------------------------------------------------------
!> @brief Perform module initialization.
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
module function plastic_phenopowerlaw_init() result(myPlasticity)
2020-08-15 19:32:10 +05:30
logical, dimension(:), allocatable :: myPlasticity
integer :: &
2021-02-16 20:36:09 +05:30
ph, i, &
2021-03-05 01:46:36 +05:30
Nmembers, &
sizeState, sizeDotState, &
startIndex, endIndex
2020-03-16 19:32:40 +05:30
integer, dimension(:), allocatable :: &
N_sl, & !< number of slip-systems for a given slip family
N_tw !< number of twin-systems for a given twin family
real(pREAL), dimension(:), allocatable :: &
2020-09-22 19:34:14 +05:30
xi_0_sl, & !< initial critical shear stress for slip
xi_0_tw, & !< initial critical shear stress for twin
2020-03-16 16:23:36 +05:30
a !< non-Schmid coefficients
character(len=:), allocatable :: &
refs, &
extmsg
type(tDict), pointer :: &
2020-08-15 19:32:10 +05:30
phases, &
phase, &
2020-11-03 03:16:46 +05:30
mech, &
2020-08-15 19:32:10 +05:30
pl
2020-08-15 19:32:10 +05:30
myPlasticity = plastic_active('phenopowerlaw')
2022-12-07 22:59:03 +05:30
if (count(myPlasticity) == 0) return
print'(/,1x,a)', '<<<+- phase:mechanical:plastic:phenopowerlaw init -+>>>'
print'(/,a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
2021-02-13 17:37:35 +05:30
phases => config_material%get_dict('phase')
2021-02-14 05:20:42 +05:30
allocate(param(phases%length))
allocate(indexDotState(phases%length))
2021-02-14 05:20:42 +05:30
allocate(state(phases%length))
extmsg = ''
2021-02-14 05:20:42 +05:30
2021-02-16 20:36:09 +05:30
do ph = 1, phases%length
if (.not. myPlasticity(ph)) cycle
2021-02-16 20:36:09 +05:30
associate(prm => param(ph), stt => state(ph), &
idx_dot => indexDotState(ph))
2021-02-16 20:36:09 +05:30
phase => phases%get_dict(ph)
mech => phase%get_dict('mechanical')
pl => mech%get_dict('plastic')
print'(/,1x,a,i0,a)', 'phase ',ph,': '//phases%key(ph)
refs = config_listReferences(pl,indent=3)
if (len(refs) > 0) print'(/,1x,a)', refs
#if defined (__GFORTRAN__)
2023-06-04 10:47:38 +05:30
prm%output = output_as1dStr(pl)
#else
2023-06-04 10:47:38 +05:30
prm%output = pl%get_as1dStr('output',defaultVal=emptyStrArray)
#endif
2018-10-14 12:56:32 +05:30
!--------------------------------------------------------------------------------------------------
! slip related parameters
2021-03-11 22:30:07 +05:30
N_sl = pl%get_as1dInt('N_sl',defaultVal=emptyIntArray)
2020-03-16 19:32:40 +05:30
prm%sum_N_sl = sum(abs(N_sl))
2020-03-16 01:46:28 +05:30
slipActive: if (prm%sum_N_sl > 0) then
prm%systems_sl = lattice_labels_slip(N_sl,phase_lattice(ph))
2021-07-21 19:16:38 +05:30
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase_lattice(ph),phase_cOverA(ph))
2021-07-21 19:16:38 +05:30
if (phase_lattice(ph) == 'cI') then
a = pl%get_as1dReal('a_nonSchmid',defaultVal=emptyRealArray)
2022-12-07 22:59:03 +05:30
if (size(a) > 0) prm%nonSchmidActive = .true.
2021-07-21 19:59:57 +05:30
prm%P_nS_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
prm%P_nS_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
else
2021-07-21 19:59:57 +05:30
prm%P_nS_pos = prm%P_sl
prm%P_nS_neg = prm%P_sl
end if
prm%h_sl_sl = lattice_interaction_SlipBySlip(N_sl,pl%get_as1dReal('h_sl-sl'),phase_lattice(ph))
xi_0_sl = pl%get_as1dReal('xi_0_sl', requiredSize=size(N_sl))
prm%xi_inf_sl = pl%get_as1dReal('xi_inf_sl', requiredSize=size(N_sl))
prm%h_int = pl%get_as1dReal('h_int', requiredSize=size(N_sl), &
defaultVal=[(0.0_pREAL,i=1,size(N_sl))])
prm%dot_gamma_0_sl = pl%get_asReal('dot_gamma_0_sl')
prm%n_sl = pl%get_asReal('n_sl')
prm%a_sl = pl%get_asReal('a_sl')
prm%h_0_sl_sl = pl%get_asReal('h_0_sl-sl')
! expand: family => system
2020-09-22 19:34:14 +05:30
xi_0_sl = math_expand(xi_0_sl, N_sl)
prm%xi_inf_sl = math_expand(prm%xi_inf_sl,N_sl)
prm%h_int = math_expand(prm%h_int, N_sl)
! sanity checks
if ( prm%dot_gamma_0_sl <= 0.0_pREAL) extmsg = trim(extmsg)//' dot_gamma_0_sl'
if ( prm%a_sl <= 0.0_pREAL) extmsg = trim(extmsg)//' a_sl'
if ( prm%n_sl <= 0.0_pREAL) extmsg = trim(extmsg)//' n_sl'
if (any(xi_0_sl <= 0.0_pREAL)) extmsg = trim(extmsg)//' xi_0_sl'
if (any(prm%xi_inf_sl <= 0.0_pREAL)) extmsg = trim(extmsg)//' xi_inf_sl'
2020-03-16 16:23:36 +05:30
else slipActive
2020-09-22 19:34:14 +05:30
xi_0_sl = emptyRealArray
allocate(prm%xi_inf_sl,prm%h_int,source=emptyRealArray)
allocate(prm%h_sl_sl(0,0))
end if slipActive
2018-10-14 12:56:32 +05:30
!--------------------------------------------------------------------------------------------------
! twin related parameters
2021-03-11 22:30:07 +05:30
N_tw = pl%get_as1dInt('N_tw', defaultVal=emptyIntArray)
2020-03-16 19:32:40 +05:30
prm%sum_N_tw = sum(abs(N_tw))
2020-03-16 01:46:28 +05:30
twinActive: if (prm%sum_N_tw > 0) then
prm%systems_tw = lattice_labels_twin(N_tw,phase_lattice(ph))
prm%P_tw = lattice_SchmidMatrix_twin(N_tw,phase_lattice(ph),phase_cOverA(ph))
prm%h_tw_tw = lattice_interaction_TwinByTwin(N_tw,pl%get_as1dReal('h_tw-tw'),phase_lattice(ph))
prm%gamma_char = lattice_characteristicShear_twin(N_tw,phase_lattice(ph),phase_cOverA(ph))
2020-09-22 19:34:14 +05:30
xi_0_tw = pl%get_as1dReal('xi_0_tw',requiredSize=size(N_tw))
2020-09-22 19:34:14 +05:30
prm%c_1 = pl%get_asReal('c_1',defaultVal=0.0_pREAL)
prm%c_2 = pl%get_asReal('c_2',defaultVal=1.0_pREAL)
prm%c_3 = pl%get_asReal('c_3',defaultVal=0.0_pREAL)
prm%c_4 = pl%get_asReal('c_4',defaultVal=0.0_pREAL)
prm%dot_gamma_0_tw = pl%get_asReal('dot_gamma_0_tw')
prm%n_tw = pl%get_asReal('n_tw')
prm%f_sat_sl_tw = pl%get_asReal('f_sat_sl-tw')
prm%h_0_tw_tw = pl%get_asReal('h_0_tw-tw')
! expand: family => system
2020-09-22 19:34:14 +05:30
xi_0_tw = math_expand(xi_0_tw,N_tw)
! sanity checks
if (prm%dot_gamma_0_tw <= 0.0_pREAL) extmsg = trim(extmsg)//' dot_gamma_0_tw'
if (prm%n_tw <= 0.0_pREAL) extmsg = trim(extmsg)//' n_tw'
2020-03-16 16:23:36 +05:30
else twinActive
2020-09-22 19:34:14 +05:30
xi_0_tw = emptyRealArray
allocate(prm%gamma_char,source=emptyRealArray)
2020-09-22 19:34:14 +05:30
allocate(prm%h_tw_tw(0,0))
end if twinActive
2018-10-14 12:56:32 +05:30
!--------------------------------------------------------------------------------------------------
! slip-twin related parameters
2020-03-16 01:46:28 +05:30
slipAndTwinActive: if (prm%sum_N_sl > 0 .and. prm%sum_N_tw > 0) then
prm%h_0_tw_sl = pl%get_asReal('h_0_tw-sl')
prm%h_sl_tw = lattice_interaction_SlipByTwin(N_sl,N_tw,pl%get_as1dReal('h_sl-tw'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
prm%h_tw_sl = lattice_interaction_TwinBySlip(N_tw,N_sl,pl%get_as1dReal('h_tw-sl'), &
2021-07-21 19:16:38 +05:30
phase_lattice(ph))
else slipAndTwinActive
2020-09-22 19:34:14 +05:30
allocate(prm%h_sl_tw(prm%sum_N_sl,prm%sum_N_tw)) ! at least one dimension is 0
allocate(prm%h_tw_sl(prm%sum_N_tw,prm%sum_N_sl)) ! at least one dimension is 0
prm%h_0_tw_sl = 0.0_pREAL
end if slipAndTwinActive
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2023-01-23 13:01:59 +05:30
Nmembers = count(material_ID_phase == ph)
sizeDotState = size(['xi_sl ','gamma_sl']) * prm%sum_N_sl &
+ size(['xi_tw ','gamma_tw']) * prm%sum_N_tw
sizeState = sizeDotState
2021-03-05 01:46:36 +05:30
call phase_allocateState(plasticState(ph),Nmembers,sizeState,sizeDotState,0)
deallocate(plasticState(ph)%dotState) ! ToDo: remove dotState completely
!--------------------------------------------------------------------------------------------------
2020-03-16 03:37:23 +05:30
! state aliases and initialization
startIndex = 1
2020-03-16 01:46:28 +05:30
endIndex = prm%sum_N_sl
idx_dot%xi_sl = [startIndex,endIndex]
stt%xi_sl => plasticState(ph)%state(startIndex:endIndex,:)
2021-07-21 18:33:28 +05:30
stt%xi_sl = spread(xi_0_sl, 2, Nmembers)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_xi',defaultVal=1.0_pREAL)
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_xi'
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_tw
idx_dot%xi_tw = [startIndex,endIndex]
stt%xi_tw => plasticState(ph)%state(startIndex:endIndex,:)
2021-07-21 18:33:28 +05:30
stt%xi_tw = spread(xi_0_tw, 2, Nmembers)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_xi',defaultVal=1.0_pREAL)
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_sl
idx_dot%gamma_sl = [startIndex,endIndex]
stt%gamma_sl => plasticState(ph)%state(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_gamma',defaultVal=1.0e-6_pREAL)
if (any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_gamma'
2020-02-11 22:06:43 +05:30
startIndex = endIndex + 1
2020-03-16 01:46:28 +05:30
endIndex = endIndex + prm%sum_N_tw
idx_dot%gamma_tw = [startIndex,endIndex]
stt%gamma_tw => plasticState(ph)%state(startIndex:endIndex,:)
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asReal('atol_gamma',defaultVal=1.0e-6_pREAL)
end associate
2020-03-15 12:58:38 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg))
2020-03-15 12:58:38 +05:30
end do
2020-08-15 19:32:10 +05:30
end function plastic_phenopowerlaw_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
2023-02-11 20:00:12 +05:30
!> @details assumes that deformation by dislocation glide affects twinned and untwinned volume
! equally (Taylor assumption). Twinning happens only in untwinned volume
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
pure module subroutine phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
Lp !< plastic velocity gradient
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) :: &
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
integer :: &
i,k,l,m,n
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_sl_pos,dot_gamma_sl_neg, &
ddot_gamma_dtau_sl_pos,ddot_gamma_dtau_sl_neg
real(pREAL), dimension(param(ph)%sum_N_tw) :: &
2021-07-23 02:52:45 +05:30
dot_gamma_tw,ddot_gamma_dtau_tw
Lp = 0.0_pREAL
dLp_dMp = 0.0_pREAL
2021-02-14 05:20:42 +05:30
associate(prm => param(ph))
2021-07-21 20:04:33 +05:30
call kinetics_sl(Mp,ph,en,dot_gamma_sl_pos,dot_gamma_sl_neg,ddot_gamma_dtau_sl_pos,ddot_gamma_dtau_sl_neg)
2020-03-16 01:46:28 +05:30
slipSystems: do i = 1, prm%sum_N_sl
2021-07-21 20:04:33 +05:30
Lp = Lp + (dot_gamma_sl_pos(i)+dot_gamma_sl_neg(i))*prm%P_sl(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
2021-07-21 20:04:33 +05:30
+ ddot_gamma_dtau_sl_pos(i) * prm%P_sl(k,l,i) * prm%P_nS_pos(m,n,i) &
+ ddot_gamma_dtau_sl_neg(i) * prm%P_sl(k,l,i) * prm%P_nS_neg(m,n,i)
end do slipSystems
2021-07-23 02:52:45 +05:30
call kinetics_tw(Mp,ph,en,dot_gamma_tw,ddot_gamma_dtau_tw)
2020-03-16 01:46:28 +05:30
twinSystems: do i = 1, prm%sum_N_tw
2021-07-21 20:04:33 +05:30
Lp = Lp + dot_gamma_tw(i)*prm%P_tw(1:3,1:3,i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = dLp_dMp(k,l,m,n) &
2021-07-23 02:52:45 +05:30
+ ddot_gamma_dtau_tw(i)*prm%P_tw(k,l,i)*prm%P_tw(m,n,i)
end do twinSystems
end associate
2021-01-26 12:03:04 +05:30
end subroutine phenopowerlaw_LpAndItsTangent
2018-08-03 04:44:18 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Calculate the rate of change of microstructure.
!--------------------------------------------------------------------------------------------------
module function phenopowerlaw_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(plasticState(ph)%sizeDotState) :: &
dotState
real(pREAL) :: &
2021-07-21 18:33:28 +05:30
xi_sl_sat_offset,&
sumF
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_sl_pos,dot_gamma_sl_neg, &
left_SlipSlip
associate(prm => param(ph), stt => state(ph), &
dot_xi_sl => dotState(indexDotState(ph)%xi_sl(1):indexDotState(ph)%xi_sl(2)), &
dot_xi_tw => dotState(indexDotState(ph)%xi_tw(1):indexDotState(ph)%xi_tw(2)), &
dot_gamma_sl => dotState(indexDotState(ph)%gamma_sl(1):indexDotState(ph)%gamma_sl(2)), &
dot_gamma_tw => dotState(indexDotState(ph)%gamma_tw(1):indexDotState(ph)%gamma_tw(2)))
2015-10-14 00:22:01 +05:30
call kinetics_sl(Mp,ph,en, dot_gamma_sl_pos,dot_gamma_sl_neg)
dot_gamma_sl = abs(dot_gamma_sl_pos+dot_gamma_sl_neg)
call kinetics_tw(Mp,ph,en, dot_gamma_tw)
2021-07-21 18:33:28 +05:30
sumF = sum(stt%gamma_tw(:,en)/prm%gamma_char)
2021-07-21 18:33:28 +05:30
xi_sl_sat_offset = prm%f_sat_sl_tw*sqrt(sumF)
left_SlipSlip = sign(abs(1.0_pREAL-stt%xi_sl(:,en) / (prm%xi_inf_sl+xi_sl_sat_offset))**prm%a_sl, &
1.0_pREAL-stt%xi_sl(:,en) / (prm%xi_inf_sl+xi_sl_sat_offset))
dot_xi_sl = prm%h_0_sl_sl * (1.0_pREAL + prm%c_1 * sumF**prm%c_2) * (1.0_pREAL + prm%h_int) &
* left_SlipSlip * matmul(prm%h_sl_sl,dot_gamma_sl) &
+ matmul(prm%h_sl_tw,dot_gamma_tw)
2021-07-21 18:33:28 +05:30
dot_xi_tw = prm%h_0_tw_sl * sum(stt%gamma_sl(:,en))**prm%c_3 &
* matmul(prm%h_tw_sl,dot_gamma_sl) &
+ prm%h_0_tw_tw * sumF**prm%c_4 * matmul(prm%h_tw_tw,dot_gamma_tw)
end associate
2018-08-26 00:43:32 +05:30
end function phenopowerlaw_dotState
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
module subroutine plastic_phenopowerlaw_result(ph,group)
2019-04-04 11:22:36 +05:30
2021-02-14 05:20:42 +05:30
integer, intent(in) :: ph
2019-04-04 11:22:36 +05:30
character(len=*), intent(in) :: group
2021-07-25 13:41:19 +05:30
integer :: ou
2019-04-04 11:22:36 +05:30
2021-02-14 05:20:42 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-25 13:41:19 +05:30
do ou = 1,size(prm%output)
select case(trim(prm%output(ou)))
case('xi_sl')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%xi_sl,group,trim(prm%output(ou)), &
'resistance against plastic slip','Pa',prm%systems_sl)
2021-07-25 13:41:19 +05:30
case('gamma_sl')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%gamma_sl,group,trim(prm%output(ou)), &
'plastic shear','1',prm%systems_sl)
2021-07-25 13:41:19 +05:30
case('xi_tw')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%xi_tw,group,trim(prm%output(ou)), &
'resistance against twinning','Pa',prm%systems_tw)
2021-07-25 13:41:19 +05:30
case('gamma_tw')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%gamma_tw,group,trim(prm%output(ou)), &
'twinning shear','1',prm%systems_tw)
2021-07-25 13:41:19 +05:30
end select
end do
2021-07-25 13:41:19 +05:30
2019-04-04 11:22:36 +05:30
end associate
2023-01-19 22:07:45 +05:30
end subroutine plastic_phenopowerlaw_result
!--------------------------------------------------------------------------------------------------
!> @brief Calculate shear rates on slip systems and their derivatives with respect to resolved
! stress.
2019-01-07 11:54:02 +05:30
!> @details Derivatives are calculated only optionally.
2023-02-11 20:00:12 +05:30
! NOTE: Contrary to common convention, here the result (i.e. intent(out)) variables have to be put
! at the end since some of them are optional.
!--------------------------------------------------------------------------------------------------
2021-07-21 18:33:28 +05:30
pure subroutine kinetics_sl(Mp,ph,en, &
2021-07-21 20:04:33 +05:30
dot_gamma_sl_pos,dot_gamma_sl_neg,ddot_gamma_dtau_sl_pos,ddot_gamma_dtau_sl_neg)
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), intent(out), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_sl_pos, &
dot_gamma_sl_neg
real(pREAL), intent(out), optional, dimension(param(ph)%sum_N_sl) :: &
2021-07-21 20:04:33 +05:30
ddot_gamma_dtau_sl_pos, &
ddot_gamma_dtau_sl_neg
real(pREAL), dimension(param(ph)%sum_N_sl) :: &
2021-07-21 18:33:28 +05:30
tau_sl_pos, &
tau_sl_neg
integer :: i
2021-02-14 05:20:42 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-23 03:39:51 +05:30
do i = 1, prm%sum_N_sl
tau_sl_pos(i) = math_tensordot(Mp,prm%P_nS_pos(1:3,1:3,i))
tau_sl_neg(i) = merge(math_tensordot(Mp,prm%P_nS_neg(1:3,1:3,i)), &
0.0_pREAL, prm%nonSchmidActive)
end do
2021-07-23 03:39:51 +05:30
where(dNeq0(tau_sl_pos))
dot_gamma_sl_pos = prm%dot_gamma_0_sl * merge(0.5_pREAL,1.0_pREAL, prm%nonSchmidActive) & ! 1/2 if non-Schmid active
2021-07-23 03:39:51 +05:30
* sign(abs(tau_sl_pos/stt%xi_sl(:,en))**prm%n_sl, tau_sl_pos)
else where
dot_gamma_sl_pos = 0.0_pREAL
end where
2021-07-23 03:39:51 +05:30
where(dNeq0(tau_sl_neg))
dot_gamma_sl_neg = prm%dot_gamma_0_sl * 0.5_pREAL & ! only used if non-Schmid active, always 1/2
2021-07-23 03:39:51 +05:30
* sign(abs(tau_sl_neg/stt%xi_sl(:,en))**prm%n_sl, tau_sl_neg)
else where
dot_gamma_sl_neg = 0.0_pREAL
end where
2021-07-23 03:39:51 +05:30
if (present(ddot_gamma_dtau_sl_pos)) then
where(dNeq0(dot_gamma_sl_pos))
ddot_gamma_dtau_sl_pos = dot_gamma_sl_pos*prm%n_sl/tau_sl_pos
else where
ddot_gamma_dtau_sl_pos = 0.0_pREAL
2021-07-23 03:39:51 +05:30
end where
end if
2021-07-23 03:39:51 +05:30
if (present(ddot_gamma_dtau_sl_neg)) then
where(dNeq0(dot_gamma_sl_neg))
ddot_gamma_dtau_sl_neg = dot_gamma_sl_neg*prm%n_sl/tau_sl_neg
else where
ddot_gamma_dtau_sl_neg = 0.0_pREAL
2021-07-23 03:39:51 +05:30
end where
end if
2021-07-23 03:39:51 +05:30
end associate
2021-07-21 18:33:28 +05:30
end subroutine kinetics_sl
!--------------------------------------------------------------------------------------------------
!> @brief Calculate shear rates on twin systems and their derivatives with respect to resolved
2023-02-11 20:00:12 +05:30
! stress. Twinning is assumed to take place only in an untwinned volume.
!> @details Derivatives are calculated and returned if corresponding output variables are present in the argument list.
! NOTE: Contrary to common convention, here the result (i.e. intent(out)) variables have to be put
! at the end since some of them are optional.
!--------------------------------------------------------------------------------------------------
2021-07-21 18:33:28 +05:30
pure subroutine kinetics_tw(Mp,ph,en,&
2021-07-21 20:04:33 +05:30
dot_gamma_tw,ddot_gamma_dtau_tw)
real(pREAL), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
2021-02-14 05:20:42 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
real(pREAL), dimension(param(ph)%sum_N_tw), intent(out) :: &
2021-07-21 20:04:33 +05:30
dot_gamma_tw
real(pREAL), dimension(param(ph)%sum_N_tw), intent(out), optional :: &
2021-07-21 20:04:33 +05:30
ddot_gamma_dtau_tw
real(pREAL), dimension(param(ph)%sum_N_tw) :: &
2021-07-21 18:33:28 +05:30
tau_tw
integer :: i
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-23 03:39:51 +05:30
tau_tw = [(math_tensordot(Mp,prm%P_tw(1:3,1:3,i)),i=1,prm%sum_N_tw)]
where(tau_tw > 0.0_pREAL)
dot_gamma_tw = (1.0_pREAL-sum(stt%gamma_tw(:,en)/prm%gamma_char)) & ! only twin in untwinned volume fraction
2021-07-23 03:39:51 +05:30
* prm%dot_gamma_0_tw*(abs(tau_tw)/stt%xi_tw(:,en))**prm%n_tw
else where
dot_gamma_tw = 0.0_pREAL
end where
2021-07-23 03:39:51 +05:30
if (present(ddot_gamma_dtau_tw)) then
where(dNeq0(dot_gamma_tw))
ddot_gamma_dtau_tw = dot_gamma_tw*prm%n_tw/tau_tw
else where
ddot_gamma_dtau_tw = 0.0_pREAL
2021-07-23 03:39:51 +05:30
end where
end if
end associate
2021-07-21 18:33:28 +05:30
end subroutine kinetics_tw
2021-01-26 05:41:32 +05:30
end submodule phenopowerlaw