2018-06-30 17:07:13 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-10-11 20:19:12 +05:30
|
|
|
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
2018-06-22 02:44:30 +05:30
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
2019-01-06 04:03:18 +05:30
|
|
|
!> @brief phenomenological crystal plasticity formulation using a powerlaw fitting
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-01-27 01:22:48 +05:30
|
|
|
submodule(phase:plastic) phenopowerlaw
|
2014-03-09 02:20:31 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
type :: tParameters
|
|
|
|
real(pReal) :: &
|
2020-09-22 22:56:39 +05:30
|
|
|
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
|
2021-03-29 10:55:55 +05:30
|
|
|
f_sat_sl_tw = 1.0_pReal, & !< push-up factor for slip saturation due to twinning
|
2020-09-22 22:56:39 +05:30
|
|
|
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
|
2020-02-14 11:47:30 +05:30
|
|
|
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)
|
2021-03-29 10:46:15 +05:30
|
|
|
gamma_char !< characteristic shear for twins
|
2020-02-14 11:47:30 +05:30
|
|
|
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
|
2020-02-14 11:47:30 +05:30
|
|
|
real(pReal), allocatable, dimension(:,:,:) :: &
|
2020-03-16 01:46:28 +05:30
|
|
|
P_sl, &
|
|
|
|
P_tw, &
|
2020-02-01 02:07:18 +05:30
|
|
|
nonSchmid_pos, &
|
|
|
|
nonSchmid_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.
|
2020-02-14 11:47:30 +05:30
|
|
|
character(len=pStringLen), allocatable, dimension(:) :: &
|
|
|
|
output
|
2020-02-01 02:07:18 +05:30
|
|
|
end type tParameters
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
type :: tPhenopowerlawState
|
|
|
|
real(pReal), pointer, dimension(:,:) :: &
|
|
|
|
xi_slip, &
|
|
|
|
xi_twin, &
|
|
|
|
gamma_slip, &
|
|
|
|
gamma_twin
|
|
|
|
end type tPhenopowerlawState
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2019-01-13 21:33:49 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! containers for parameters and state
|
2020-02-01 02:07:18 +05:30
|
|
|
type(tParameters), allocatable, dimension(:) :: param
|
|
|
|
type(tPhenopowerlawState), allocatable, dimension(:) :: &
|
|
|
|
dotState, &
|
|
|
|
state
|
2015-10-30 21:18:30 +05:30
|
|
|
|
2012-03-09 01:55:28 +05:30
|
|
|
contains
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2013-07-01 11:40:42 +05:30
|
|
|
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 11:47:30 +05:30
|
|
|
!> @brief Perform module initialization.
|
2013-07-01 11:40:42 +05:30
|
|
|
!> @details reads in material parameters, allocates arrays, and does sanity checks
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-09-14 00:30:34 +05:30
|
|
|
module function plastic_phenopowerlaw_init() result(myPlasticity)
|
2009-10-21 18:40:12 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
logical, dimension(:), allocatable :: myPlasticity
|
2020-02-01 02:07:18 +05:30
|
|
|
integer :: &
|
2021-02-16 20:36:09 +05:30
|
|
|
ph, i, &
|
2021-03-05 01:46:36 +05:30
|
|
|
Nmembers, &
|
2020-02-01 02:07:18 +05:30
|
|
|
sizeState, sizeDotState, &
|
|
|
|
startIndex, endIndex
|
2020-03-16 19:32:40 +05:30
|
|
|
integer, dimension(:), allocatable :: &
|
2020-03-16 03:37:23 +05:30
|
|
|
N_sl, N_tw
|
2020-03-16 19:32:40 +05:30
|
|
|
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
|
2020-02-01 02:07:18 +05:30
|
|
|
character(len=pStringLen) :: &
|
|
|
|
extmsg = ''
|
2020-08-15 19:32:10 +05:30
|
|
|
class(tNode), pointer :: &
|
|
|
|
phases, &
|
|
|
|
phase, &
|
2020-11-03 03:16:46 +05:30
|
|
|
mech, &
|
2020-08-15 19:32:10 +05:30
|
|
|
pl
|
2020-02-14 11:47:30 +05:30
|
|
|
|
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
myPlasticity = plastic_active('phenopowerlaw')
|
2021-02-14 05:20:42 +05:30
|
|
|
if(count(myPlasticity) == 0) return
|
2020-09-14 00:30:34 +05:30
|
|
|
|
2021-02-13 23:27:41 +05:30
|
|
|
print'(/,a)', ' <<<+- phase:mechanical:plastic:phenopowerlaw init -+>>>'
|
2021-02-14 05:20:42 +05:30
|
|
|
print'(a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
|
2021-02-13 17:37:35 +05:30
|
|
|
|
|
|
|
|
2020-09-13 14:09:17 +05:30
|
|
|
phases => config_material%get('phase')
|
2021-02-14 05:20:42 +05:30
|
|
|
allocate(param(phases%length))
|
|
|
|
allocate(state(phases%length))
|
|
|
|
allocate(dotState(phases%length))
|
|
|
|
|
2021-02-16 20:36:09 +05:30
|
|
|
do ph = 1, phases%length
|
|
|
|
if(.not. myPlasticity(ph)) cycle
|
|
|
|
|
|
|
|
associate(prm => param(ph), dot => dotState(ph), stt => state(ph))
|
|
|
|
|
|
|
|
phase => phases%get(ph)
|
2021-03-25 23:52:59 +05:30
|
|
|
mech => phase%get('mechanical')
|
|
|
|
pl => mech%get('plastic')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
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
|
2020-08-15 19:32:10 +05:30
|
|
|
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
|
|
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-11-29 16:28:39 +05:30
|
|
|
if(phase%get_asString('lattice') == 'cI') then
|
2021-03-11 22:30:07 +05:30
|
|
|
a = pl%get_as1dFloat('a_nonSchmid',defaultVal=emptyRealArray)
|
2020-03-16 16:23:36 +05:30
|
|
|
if(size(a) > 0) prm%nonSchmidActive = .true.
|
|
|
|
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
|
|
|
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
|
2020-02-01 02:07:18 +05:30
|
|
|
else
|
2020-03-16 01:46:28 +05:30
|
|
|
prm%nonSchmid_pos = prm%P_sl
|
|
|
|
prm%nonSchmid_neg = prm%P_sl
|
2020-02-01 02:07:18 +05:30
|
|
|
endif
|
2020-09-22 19:34:14 +05:30
|
|
|
prm%h_sl_sl = lattice_interaction_SlipBySlip(N_sl, &
|
2021-06-19 18:43:27 +05:30
|
|
|
pl%get_as1dFloat('h_sl-sl'), &
|
2020-09-22 19:34:14 +05:30
|
|
|
phase%get_asString('lattice'))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-03-11 22:30:07 +05:30
|
|
|
xi_0_sl = pl%get_as1dFloat('xi_0_sl', requiredSize=size(N_sl))
|
|
|
|
prm%xi_inf_sl = pl%get_as1dFloat('xi_inf_sl', requiredSize=size(N_sl))
|
|
|
|
prm%h_int = pl%get_as1dFloat('h_int', requiredSize=size(N_sl), &
|
2021-02-16 20:36:09 +05:30
|
|
|
defaultVal=[(0.0_pReal,i=1,size(N_sl))])
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-09-22 19:34:14 +05:30
|
|
|
prm%dot_gamma_0_sl = pl%get_asFloat('dot_gamma_0_sl')
|
|
|
|
prm%n_sl = pl%get_asFloat('n_sl')
|
|
|
|
prm%a_sl = pl%get_asFloat('a_sl')
|
2021-06-19 18:43:27 +05:30
|
|
|
prm%h_0_sl_sl = pl%get_asFloat('h_0_sl-sl')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! 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)
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! sanity checks
|
2020-09-22 19:34:14 +05:30
|
|
|
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
|
|
|
|
2020-02-01 02:07:18 +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))
|
2020-02-01 02:07:18 +05:30
|
|
|
endif slipActive
|
2020-02-14 11:47:30 +05:30
|
|
|
|
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
|
2020-09-22 19:34:14 +05:30
|
|
|
prm%P_tw = lattice_SchmidMatrix_twin(N_tw,phase%get_asString('lattice'),&
|
|
|
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
|
|
|
prm%h_tw_tw = lattice_interaction_TwinByTwin(N_tw,&
|
2021-06-19 18:43:27 +05:30
|
|
|
pl%get_as1dFloat('h_tw-tw'), &
|
2020-09-22 19:34:14 +05:30
|
|
|
phase%get_asString('lattice'))
|
2021-03-29 10:46:15 +05:30
|
|
|
prm%gamma_char = lattice_characteristicShear_twin(N_tw,phase%get_asString('lattice'),&
|
2020-09-22 19:34:14 +05:30
|
|
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
|
|
|
|
2021-03-11 22:30:07 +05:30
|
|
|
xi_0_tw = pl%get_as1dFloat('xi_0_tw',requiredSize=size(N_tw))
|
2020-09-22 19:34:14 +05:30
|
|
|
|
|
|
|
prm%c_1 = pl%get_asFloat('c_1',defaultVal=0.0_pReal)
|
|
|
|
prm%c_2 = pl%get_asFloat('c_2',defaultVal=1.0_pReal)
|
|
|
|
prm%c_3 = pl%get_asFloat('c_3',defaultVal=0.0_pReal)
|
|
|
|
prm%c_4 = pl%get_asFloat('c_4',defaultVal=0.0_pReal)
|
|
|
|
prm%dot_gamma_0_tw = pl%get_asFloat('dot_gamma_0_tw')
|
|
|
|
prm%n_tw = pl%get_asFloat('n_tw')
|
2021-06-19 18:43:27 +05:30
|
|
|
prm%f_sat_sl_tw = pl%get_asFloat('f_sat_sl-tw')
|
|
|
|
prm%h_0_tw_tw = pl%get_asFloat('h_0_tw-tw')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! expand: family => system
|
2020-09-22 19:34:14 +05:30
|
|
|
xi_0_tw = math_expand(xi_0_tw,N_tw)
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! sanity checks
|
2020-09-22 19:34:14 +05:30
|
|
|
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
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
else twinActive
|
2020-09-22 19:34:14 +05:30
|
|
|
xi_0_tw = emptyRealArray
|
2021-03-29 10:46:15 +05:30
|
|
|
allocate(prm%gamma_char,source=emptyRealArray)
|
2020-09-22 19:34:14 +05:30
|
|
|
allocate(prm%h_tw_tw(0,0))
|
2020-02-01 02:07:18 +05:30
|
|
|
endif twinActive
|
2020-02-14 11:47:30 +05:30
|
|
|
|
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
|
2021-06-19 18:43:27 +05:30
|
|
|
prm%h_0_tw_sl = pl%get_asFloat('h_0_tw-sl')
|
2020-09-22 19:34:14 +05:30
|
|
|
prm%h_sl_tw = lattice_interaction_SlipByTwin(N_sl,N_tw,&
|
2021-06-19 18:43:27 +05:30
|
|
|
pl%get_as1dFloat('h_sl-tw'), &
|
2020-09-22 19:34:14 +05:30
|
|
|
phase%get_asString('lattice'))
|
|
|
|
prm%h_tw_sl = lattice_interaction_TwinBySlip(N_tw,N_sl,&
|
2021-06-19 18:43:27 +05:30
|
|
|
pl%get_as1dFloat('h_tw-sl'), &
|
2020-09-22 19:34:14 +05:30
|
|
|
phase%get_asString('lattice'))
|
2020-02-01 02:07:18 +05:30
|
|
|
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
|
2020-02-01 02:07:18 +05:30
|
|
|
endif slipAndTwinActive
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-10-14 12:56:32 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! output pararameters
|
2020-08-15 19:32:10 +05:30
|
|
|
|
|
|
|
#if defined (__GFORTRAN__)
|
2021-03-11 22:30:07 +05:30
|
|
|
prm%output = output_as1dString(pl)
|
2020-08-15 19:32:10 +05:30
|
|
|
#else
|
2021-03-11 22:30:07 +05:30
|
|
|
prm%output = pl%get_as1dString('output',defaultVal=emptyStringArray)
|
2020-08-15 19:32:10 +05:30
|
|
|
#endif
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2014-07-02 17:57:39 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! allocate state arrays
|
2021-04-06 15:08:44 +05:30
|
|
|
Nmembers = count(material_phaseID == ph)
|
2020-03-16 04:57:52 +05:30
|
|
|
sizeDotState = size(['xi_sl ','gamma_sl']) * prm%sum_N_sl &
|
|
|
|
+ size(['xi_tw ','gamma_tw']) * prm%sum_N_tw
|
2020-02-01 02:07:18 +05:30
|
|
|
sizeState = sizeDotState
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-09-14 00:30:34 +05:30
|
|
|
|
2021-03-05 01:46:36 +05:30
|
|
|
call phase_allocateState(plasticState(ph),Nmembers,sizeState,sizeDotState,0)
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-04-25 23:11:18 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-03-16 03:37:23 +05:30
|
|
|
! state aliases and initialization
|
2020-02-01 02:07:18 +05:30
|
|
|
startIndex = 1
|
2020-03-16 01:46:28 +05:30
|
|
|
endIndex = prm%sum_N_sl
|
2021-02-16 20:36:09 +05:30
|
|
|
stt%xi_slip => plasticState(ph)%state (startIndex:endIndex,:)
|
2021-03-05 01:46:36 +05:30
|
|
|
stt%xi_slip = spread(xi_0_sl, 2, Nmembers)
|
2021-02-16 20:36:09 +05:30
|
|
|
dot%xi_slip => plasticState(ph)%dotState(startIndex:endIndex,:)
|
|
|
|
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_xi',defaultVal=1.0_pReal)
|
|
|
|
if(any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_xi'
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
startIndex = endIndex + 1
|
2020-03-16 01:46:28 +05:30
|
|
|
endIndex = endIndex + prm%sum_N_tw
|
2021-02-16 20:36:09 +05:30
|
|
|
stt%xi_twin => plasticState(ph)%state (startIndex:endIndex,:)
|
2021-03-05 01:46:36 +05:30
|
|
|
stt%xi_twin = spread(xi_0_tw, 2, Nmembers)
|
2021-02-16 20:36:09 +05:30
|
|
|
dot%xi_twin => plasticState(ph)%dotState(startIndex:endIndex,:)
|
|
|
|
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_xi',defaultVal=1.0_pReal)
|
|
|
|
if(any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_xi'
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
startIndex = endIndex + 1
|
2020-03-16 01:46:28 +05:30
|
|
|
endIndex = endIndex + prm%sum_N_sl
|
2021-02-16 20:36:09 +05:30
|
|
|
stt%gamma_slip => plasticState(ph)%state (startIndex:endIndex,:)
|
|
|
|
dot%gamma_slip => plasticState(ph)%dotState(startIndex:endIndex,:)
|
|
|
|
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_gamma',defaultVal=1.0e-6_pReal)
|
|
|
|
if(any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_gamma'
|
2020-02-01 02:07:18 +05:30
|
|
|
! global alias
|
2021-02-16 20:36:09 +05:30
|
|
|
plasticState(ph)%slipRate => plasticState(ph)%dotState(startIndex:endIndex,:)
|
2020-02-11 22:06:43 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
startIndex = endIndex + 1
|
2020-03-16 01:46:28 +05:30
|
|
|
endIndex = endIndex + prm%sum_N_tw
|
2021-02-16 20:36:09 +05:30
|
|
|
stt%gamma_twin => plasticState(ph)%state (startIndex:endIndex,:)
|
|
|
|
dot%gamma_twin => plasticState(ph)%dotState(startIndex:endIndex,:)
|
|
|
|
plasticState(ph)%atol(startIndex:endIndex) = pl%get_asFloat('atol_gamma',defaultVal=1.0e-6_pReal)
|
|
|
|
if(any(plasticState(ph)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_gamma'
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
end associate
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-15 12:58:38 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! exit if any parameter is out of range
|
2020-08-15 19:32:10 +05:30
|
|
|
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(phenopowerlaw)')
|
2020-03-15 12:58:38 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
end function plastic_phenopowerlaw_init
|
2009-07-22 21:37:19 +05:30
|
|
|
|
|
|
|
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 11:47:30 +05:30
|
|
|
!> @brief Calculate plastic velocity gradient and its tangent.
|
2019-01-07 11:54:02 +05:30
|
|
|
!> @details asummes that deformation by dislocation glide affects twinned and untwinned volume
|
2019-01-07 01:48:35 +05:30
|
|
|
! equally (Taylor assumption). Twinning happens only in untwinned volume
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-04-25 11:36:52 +05:30
|
|
|
pure module subroutine phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
|
2018-10-14 19:23:24 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
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
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
real(pReal), dimension(3,3), intent(in) :: &
|
|
|
|
Mp !< Mandel stress
|
|
|
|
integer, intent(in) :: &
|
2021-01-27 03:11:41 +05:30
|
|
|
ph, &
|
2021-04-25 11:36:52 +05:30
|
|
|
en
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
integer :: &
|
|
|
|
i,k,l,m,n
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_slip_pos,gdot_slip_neg, &
|
|
|
|
dgdot_dtauslip_pos,dgdot_dtauslip_neg
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_tw) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_twin,dgdot_dtautwin
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
Lp = 0.0_pReal
|
|
|
|
dLp_dMp = 0.0_pReal
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
associate(prm => param(ph))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-04-25 11:36:52 +05:30
|
|
|
call kinetics_slip(Mp,ph,en,gdot_slip_pos,gdot_slip_neg,dgdot_dtauslip_pos,dgdot_dtauslip_neg)
|
2020-03-16 01:46:28 +05:30
|
|
|
slipSystems: do i = 1, prm%sum_N_sl
|
|
|
|
Lp = Lp + (gdot_slip_pos(i)+gdot_slip_neg(i))*prm%P_sl(1:3,1:3,i)
|
2020-02-01 02:07:18 +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) &
|
2020-03-16 01:46:28 +05:30
|
|
|
+ dgdot_dtauslip_pos(i) * prm%P_sl(k,l,i) * prm%nonSchmid_pos(m,n,i) &
|
|
|
|
+ dgdot_dtauslip_neg(i) * prm%P_sl(k,l,i) * prm%nonSchmid_neg(m,n,i)
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo slipSystems
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-04-25 11:36:52 +05:30
|
|
|
call kinetics_twin(Mp,ph,en,gdot_twin,dgdot_dtautwin)
|
2020-03-16 01:46:28 +05:30
|
|
|
twinSystems: do i = 1, prm%sum_N_tw
|
|
|
|
Lp = Lp + gdot_twin(i)*prm%P_tw(1:3,1:3,i)
|
2020-02-01 02:07:18 +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) &
|
2020-03-16 01:46:28 +05:30
|
|
|
+ dgdot_dtautwin(i)*prm%P_tw(k,l,i)*prm%P_tw(m,n,i)
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo twinSystems
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
end associate
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2021-01-26 12:03:04 +05:30
|
|
|
end subroutine phenopowerlaw_LpAndItsTangent
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2018-08-03 04:44:18 +05:30
|
|
|
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 11:47:30 +05:30
|
|
|
!> @brief Calculate the rate of change of microstructure.
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-04-25 11:36:52 +05:30
|
|
|
module subroutine phenopowerlaw_dotState(Mp,ph,en)
|
2018-10-14 19:23:24 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
real(pReal), dimension(3,3), intent(in) :: &
|
|
|
|
Mp !< Mandel stress
|
|
|
|
integer, intent(in) :: &
|
2021-01-27 03:11:41 +05:30
|
|
|
ph, &
|
2021-04-25 11:36:52 +05:30
|
|
|
en
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
real(pReal) :: &
|
|
|
|
c_SlipSlip,c_TwinSlip,c_TwinTwin, &
|
|
|
|
xi_slip_sat_offset,&
|
|
|
|
sumGamma,sumF
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
left_SlipSlip,right_SlipSlip, &
|
|
|
|
gdot_slip_pos,gdot_slip_neg
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
associate(prm => param(ph), stt => state(ph), &
|
|
|
|
dot => dotState(ph))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-04-25 11:36:52 +05:30
|
|
|
sumGamma = sum(stt%gamma_slip(:,en))
|
|
|
|
sumF = sum(stt%gamma_twin(:,en)/prm%gamma_char)
|
2015-10-14 00:22:01 +05:30
|
|
|
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-10-22 20:25:07 +05:30
|
|
|
! system-independent (nonlinear) prefactors to M_Xx (X influenced by x) matrices
|
2020-09-22 19:34:14 +05:30
|
|
|
c_SlipSlip = prm%h_0_sl_sl * (1.0_pReal + prm%c_1*sumF** prm%c_2)
|
|
|
|
c_TwinSlip = prm%h_0_tw_sl * sumGamma**prm%c_3
|
|
|
|
c_TwinTwin = prm%h_0_tw_tw * sumF**prm%c_4
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2013-07-01 11:40:42 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-26 00:25:15 +05:30
|
|
|
! calculate left and right vectors
|
2020-09-22 19:34:14 +05:30
|
|
|
left_SlipSlip = 1.0_pReal + prm%h_int
|
2021-03-29 10:55:55 +05:30
|
|
|
xi_slip_sat_offset = prm%f_sat_sl_tw*sqrt(sumF)
|
2021-04-25 11:36:52 +05:30
|
|
|
right_SlipSlip = abs(1.0_pReal-stt%xi_slip(:,en) / (prm%xi_inf_sl+xi_slip_sat_offset)) **prm%a_sl &
|
|
|
|
* sign(1.0_pReal,1.0_pReal-stt%xi_slip(:,en) / (prm%xi_inf_sl+xi_slip_sat_offset))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-09-12 13:29:09 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! shear rates
|
2021-04-25 11:36:52 +05:30
|
|
|
call kinetics_slip(Mp,ph,en,gdot_slip_pos,gdot_slip_neg)
|
|
|
|
dot%gamma_slip(:,en) = abs(gdot_slip_pos+gdot_slip_neg)
|
|
|
|
call kinetics_twin(Mp,ph,en,dot%gamma_twin(:,en))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2012-10-11 20:19:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-09-12 12:53:11 +05:30
|
|
|
! hardening
|
2021-04-25 11:36:52 +05:30
|
|
|
dot%xi_slip(:,en) = c_SlipSlip * left_SlipSlip * &
|
|
|
|
matmul(prm%h_sl_sl,dot%gamma_slip(:,en)*right_SlipSlip) &
|
|
|
|
+ matmul(prm%h_sl_tw,dot%gamma_twin(:,en))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-04-25 11:36:52 +05:30
|
|
|
dot%xi_twin(:,en) = c_TwinSlip * matmul(prm%h_tw_sl,dot%gamma_slip(:,en)) &
|
|
|
|
+ c_TwinTwin * matmul(prm%h_tw_tw,dot%gamma_twin(:,en))
|
2020-02-01 02:07:18 +05:30
|
|
|
end associate
|
2018-08-26 00:43:32 +05:30
|
|
|
|
2021-01-26 16:47:00 +05:30
|
|
|
end subroutine phenopowerlaw_dotState
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2018-06-30 17:41:05 +05:30
|
|
|
|
2019-01-07 01:26:36 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 11:47:30 +05:30
|
|
|
!> @brief Write results to HDF5 output file.
|
2019-01-07 01:26:36 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-02-14 05:20:42 +05:30
|
|
|
module subroutine plastic_phenopowerlaw_results(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
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2019-04-04 11:22:36 +05:30
|
|
|
integer :: o
|
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
associate(prm => param(ph), stt => state(ph))
|
2020-02-14 11:47:30 +05:30
|
|
|
outputsLoop: do o = 1,size(prm%output)
|
|
|
|
select case(trim(prm%output(o)))
|
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
case('xi_sl')
|
2021-06-01 20:35:13 +05:30
|
|
|
if(prm%sum_N_sl>0) call results_writeDataset(stt%xi_slip,group,trim(prm%output(o)), &
|
2020-03-16 01:46:28 +05:30
|
|
|
'resistance against plastic slip','Pa')
|
2020-08-15 19:32:10 +05:30
|
|
|
case('gamma_sl')
|
2021-06-01 20:35:13 +05:30
|
|
|
if(prm%sum_N_sl>0) call results_writeDataset(stt%gamma_slip,group,trim(prm%output(o)), &
|
2020-03-16 01:46:28 +05:30
|
|
|
'plastic shear','1')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
case('xi_tw')
|
2021-06-01 20:35:13 +05:30
|
|
|
if(prm%sum_N_tw>0) call results_writeDataset(stt%xi_twin,group,trim(prm%output(o)), &
|
2020-03-16 01:46:28 +05:30
|
|
|
'resistance against twinning','Pa')
|
2020-08-15 19:32:10 +05:30
|
|
|
case('gamma_tw')
|
2021-06-01 20:35:13 +05:30
|
|
|
if(prm%sum_N_tw>0) call results_writeDataset(stt%gamma_twin,group,trim(prm%output(o)), &
|
2020-03-16 01:46:28 +05:30
|
|
|
'twinning shear','1')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2019-04-04 11:22:36 +05:30
|
|
|
end select
|
|
|
|
enddo outputsLoop
|
|
|
|
end associate
|
|
|
|
|
2019-01-07 01:26:36 +05:30
|
|
|
end subroutine plastic_phenopowerlaw_results
|
|
|
|
|
|
|
|
|
2018-09-12 12:53:11 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 13:30:14 +05:30
|
|
|
!> @brief Calculate shear rates on slip systems and their derivatives with respect to resolved
|
2020-02-14 11:47:30 +05:30
|
|
|
! stress.
|
2019-01-07 11:54:02 +05:30
|
|
|
!> @details Derivatives are calculated only optionally.
|
|
|
|
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
|
2020-02-14 11:47:30 +05:30
|
|
|
! have the optional arguments at the end.
|
2018-09-12 12:53:11 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-04-25 11:36:52 +05:30
|
|
|
pure subroutine kinetics_slip(Mp,ph,en, &
|
2018-12-21 19:26:32 +05:30
|
|
|
gdot_slip_pos,gdot_slip_neg,dgdot_dtau_slip_pos,dgdot_dtau_slip_neg)
|
2018-09-12 12:53:11 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
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
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), intent(out), dimension(param(ph)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_slip_pos, &
|
|
|
|
gdot_slip_neg
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), intent(out), optional, dimension(param(ph)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
dgdot_dtau_slip_pos, &
|
|
|
|
dgdot_dtau_slip_neg
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
tau_slip_pos, &
|
|
|
|
tau_slip_neg
|
|
|
|
integer :: i
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
associate(prm => param(ph), stt => state(ph))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
do i = 1, prm%sum_N_sl
|
2020-03-16 23:13:04 +05:30
|
|
|
tau_slip_pos(i) = math_tensordot(Mp,prm%nonSchmid_pos(1:3,1:3,i))
|
|
|
|
tau_slip_neg(i) = merge(math_tensordot(Mp,prm%nonSchmid_neg(1:3,1:3,i)), &
|
2020-03-16 16:23:36 +05:30
|
|
|
0.0_pReal, prm%nonSchmidActive)
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
where(dNeq0(tau_slip_pos))
|
2020-09-22 19:34:14 +05:30
|
|
|
gdot_slip_pos = prm%dot_gamma_0_sl * merge(0.5_pReal,1.0_pReal, prm%nonSchmidActive) & ! 1/2 if non-Schmid active
|
2021-04-25 11:36:52 +05:30
|
|
|
* sign(abs(tau_slip_pos/stt%xi_slip(:,en))**prm%n_sl, tau_slip_pos)
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
gdot_slip_pos = 0.0_pReal
|
|
|
|
end where
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
where(dNeq0(tau_slip_neg))
|
2020-09-22 19:34:14 +05:30
|
|
|
gdot_slip_neg = prm%dot_gamma_0_sl * 0.5_pReal & ! only used if non-Schmid active, always 1/2
|
2021-04-25 11:36:52 +05:30
|
|
|
* sign(abs(tau_slip_neg/stt%xi_slip(:,en))**prm%n_sl, tau_slip_neg)
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
gdot_slip_neg = 0.0_pReal
|
|
|
|
end where
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
if (present(dgdot_dtau_slip_pos)) then
|
|
|
|
where(dNeq0(gdot_slip_pos))
|
2020-09-22 19:34:14 +05:30
|
|
|
dgdot_dtau_slip_pos = gdot_slip_pos*prm%n_sl/tau_slip_pos
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
dgdot_dtau_slip_pos = 0.0_pReal
|
|
|
|
end where
|
|
|
|
endif
|
|
|
|
if (present(dgdot_dtau_slip_neg)) then
|
|
|
|
where(dNeq0(gdot_slip_neg))
|
2020-09-22 19:34:14 +05:30
|
|
|
dgdot_dtau_slip_neg = gdot_slip_neg*prm%n_sl/tau_slip_neg
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
dgdot_dtau_slip_neg = 0.0_pReal
|
|
|
|
end where
|
|
|
|
endif
|
|
|
|
end associate
|
2018-09-12 16:35:23 +05:30
|
|
|
|
|
|
|
end subroutine kinetics_slip
|
2018-09-12 12:53:11 +05:30
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-14 11:47:30 +05:30
|
|
|
!> @brief Calculate shear rates on twin systems and their derivatives with respect to resolved
|
|
|
|
! stress. Twinning is assumed to take place only in untwinned volume.
|
|
|
|
!> @details Derivatives are calculated only optionally.
|
2019-01-07 11:54:02 +05:30
|
|
|
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
|
2018-12-12 20:43:57 +05:30
|
|
|
! have the optional arguments at the end.
|
2018-09-12 12:53:11 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-04-25 11:36:52 +05:30
|
|
|
pure subroutine kinetics_twin(Mp,ph,en,&
|
2018-12-21 19:26:32 +05:30
|
|
|
gdot_twin,dgdot_dtau_twin)
|
2018-09-12 12:53:11 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
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
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_tw), intent(out) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_twin
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_tw), intent(out), optional :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
dgdot_dtau_twin
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
real(pReal), dimension(param(ph)%sum_N_tw) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
tau_twin
|
|
|
|
integer :: i
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2021-02-14 05:20:42 +05:30
|
|
|
associate(prm => param(ph), stt => state(ph))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
do i = 1, prm%sum_N_tw
|
2020-03-16 23:13:04 +05:30
|
|
|
tau_twin(i) = math_tensordot(Mp,prm%P_tw(1:3,1:3,i))
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
where(tau_twin > 0.0_pReal)
|
2021-04-25 11:36:52 +05:30
|
|
|
gdot_twin = (1.0_pReal-sum(stt%gamma_twin(:,en)/prm%gamma_char)) & ! only twin in untwinned volume fraction
|
|
|
|
* prm%dot_gamma_0_tw*(abs(tau_twin)/stt%xi_twin(:,en))**prm%n_tw
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
gdot_twin = 0.0_pReal
|
|
|
|
end where
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
if (present(dgdot_dtau_twin)) then
|
|
|
|
where(dNeq0(gdot_twin))
|
2020-09-22 19:34:14 +05:30
|
|
|
dgdot_dtau_twin = gdot_twin*prm%n_tw/tau_twin
|
2020-02-01 02:07:18 +05:30
|
|
|
else where
|
|
|
|
dgdot_dtau_twin = 0.0_pReal
|
|
|
|
end where
|
|
|
|
endif
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
end associate
|
2018-08-26 00:57:14 +05:30
|
|
|
|
2018-12-21 19:26:32 +05:30
|
|
|
end subroutine kinetics_twin
|
2018-12-12 11:10:57 +05:30
|
|
|
|
2021-01-26 05:41:32 +05:30
|
|
|
end submodule phenopowerlaw
|