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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-03 01:13:02 +05:30
|
|
|
submodule(constitutive) plastic_phenopowerlaw
|
2014-03-09 02:20:31 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
type :: tParameters
|
|
|
|
real(pReal) :: &
|
2020-03-16 16:23:36 +05:30
|
|
|
gdot0_slip = 1.0_pReal, & !< reference shear strain rate for slip
|
|
|
|
gdot0_twin = 1.0_pReal, & !< reference shear strain rate for twin
|
|
|
|
n_slip = 1.0_pReal, & !< stress exponent for slip
|
|
|
|
n_twin = 1.0_pReal, & !< stress exponent for twin
|
|
|
|
spr = 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, &
|
|
|
|
h0_SlipSlip = 1.0_pReal, & !< reference hardening slip - slip
|
|
|
|
h0_TwinSlip = 1.0_pReal, & !< reference hardening twin - slip
|
|
|
|
h0_TwinTwin = 1.0_pReal, & !< reference hardening twin - twin
|
|
|
|
a_slip = 1.0_pReal
|
2020-02-14 11:47:30 +05:30
|
|
|
real(pReal), allocatable, dimension(:) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
xi_slip_sat, & !< maximum critical shear stress for slip
|
|
|
|
H_int, & !< per family hardening activity (optional)
|
|
|
|
gamma_twin_char !< characteristic shear for twins
|
2020-02-14 11:47:30 +05:30
|
|
|
real(pReal), allocatable, dimension(:,:) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
interaction_SlipSlip, & !< slip resistance from slip activity
|
|
|
|
interaction_SlipTwin, & !< slip resistance from twin activity
|
|
|
|
interaction_TwinSlip, & !< twin resistance from slip activity
|
|
|
|
interaction_TwinTwin !< 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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-03 01:13:02 +05:30
|
|
|
module subroutine plastic_phenopowerlaw_init
|
2009-10-21 18:40:12 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
integer :: &
|
|
|
|
Ninstance, &
|
|
|
|
p, i, &
|
2020-02-14 11:47:30 +05:30
|
|
|
NipcMyPhase, &
|
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-03-16 16:23:36 +05:30
|
|
|
xi_slip_0, & !< initial critical shear stress for slip
|
|
|
|
xi_twin_0, & !< initial critical shear stress for twin
|
|
|
|
a !< non-Schmid coefficients
|
2020-02-01 02:07:18 +05:30
|
|
|
character(len=pStringLen) :: &
|
|
|
|
extmsg = ''
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-15 12:58:38 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- plastic_'//PLASTICITY_PHENOPOWERLAW_LABEL//' init -+>>>'; flush(6)
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
Ninstance = count(phase_plasticity == PLASTICITY_PHENOPOWERLAW_ID)
|
|
|
|
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
|
|
|
|
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
allocate(param(Ninstance))
|
|
|
|
allocate(state(Ninstance))
|
|
|
|
allocate(dotState(Ninstance))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
do p = 1, size(phase_plasticity)
|
|
|
|
if (phase_plasticity(p) /= PLASTICITY_PHENOPOWERLAW_ID) cycle
|
|
|
|
associate(prm => param(phase_plasticityInstance(p)), &
|
|
|
|
dot => dotState(phase_plasticityInstance(p)), &
|
|
|
|
stt => state(phase_plasticityInstance(p)), &
|
|
|
|
config => config_phase(p))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-10-14 12:56:32 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! slip related parameters
|
2020-03-16 03:37:23 +05:30
|
|
|
N_sl = config%getInts('nslip',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-03-16 03:37:23 +05:30
|
|
|
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,config%getString('lattice_structure'),&
|
|
|
|
config%getFloat('c/a',defaultVal=0.0_pReal))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
if(trim(config%getString('lattice_structure')) == 'bcc') then
|
2020-03-16 16:23:36 +05:30
|
|
|
a = config%getFloats('nonschmid_coefficients',defaultVal = emptyRealArray)
|
|
|
|
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-03-16 03:37:23 +05:30
|
|
|
prm%interaction_SlipSlip = lattice_interaction_SlipBySlip(N_sl, &
|
2020-02-01 02:07:18 +05:30
|
|
|
config%getFloats('interaction_slipslip'), &
|
|
|
|
config%getString('lattice_structure'))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_slip_0 = config%getFloats('tau0_slip', requiredSize=size(N_sl))
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%xi_slip_sat = config%getFloats('tausat_slip', requiredSize=size(N_sl))
|
|
|
|
prm%H_int = config%getFloats('h_int', requiredSize=size(N_sl), &
|
|
|
|
defaultVal=[(0.0_pReal,i=1,size(N_sl))])
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
prm%gdot0_slip = config%getFloat('gdot0_slip')
|
|
|
|
prm%n_slip = config%getFloat('n_slip')
|
|
|
|
prm%a_slip = config%getFloat('a_slip')
|
|
|
|
prm%h0_SlipSlip = config%getFloat('h0_slipslip')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! expand: family => system
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_slip_0 = math_expand(xi_slip_0, N_sl)
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%xi_slip_sat = math_expand(prm%xi_slip_sat,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
|
|
|
|
if ( prm%gdot0_slip <= 0.0_pReal) extmsg = trim(extmsg)//' gdot0_slip'
|
|
|
|
if ( prm%a_slip <= 0.0_pReal) extmsg = trim(extmsg)//' a_slip'
|
|
|
|
if ( prm%n_slip <= 0.0_pReal) extmsg = trim(extmsg)//' n_slip'
|
2020-03-16 16:23:36 +05:30
|
|
|
if (any(xi_slip_0 <= 0.0_pReal)) extmsg = trim(extmsg)//' xi_slip_0'
|
2020-02-01 02:07:18 +05:30
|
|
|
if (any(prm%xi_slip_sat <= 0.0_pReal)) extmsg = trim(extmsg)//' xi_slip_sat'
|
2020-03-16 16:23:36 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
else slipActive
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_slip_0 = emptyRealArray
|
|
|
|
allocate(prm%xi_slip_sat,prm%H_int,source=emptyRealArray)
|
2020-02-01 02:07:18 +05:30
|
|
|
allocate(prm%interaction_SlipSlip(0,0))
|
|
|
|
endif slipActive
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-10-14 12:56:32 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! twin related parameters
|
2020-03-16 03:37:23 +05:30
|
|
|
N_tw = config%getInts('ntwin', 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-03-16 03:37:23 +05:30
|
|
|
prm%P_tw = lattice_SchmidMatrix_twin(N_tw,config%getString('lattice_structure'),&
|
2020-02-01 02:07:18 +05:30
|
|
|
config%getFloat('c/a',defaultVal=0.0_pReal))
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%interaction_TwinTwin = lattice_interaction_TwinByTwin(N_tw,&
|
2020-02-01 02:07:18 +05:30
|
|
|
config%getFloats('interaction_twintwin'), &
|
|
|
|
config%getString('lattice_structure'))
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%gamma_twin_char = lattice_characteristicShear_twin(N_tw,config%getString('lattice_structure'),&
|
2020-04-28 22:37:17 +05:30
|
|
|
config%getFloat('c/a',defaultVal=0.0_pReal))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_twin_0 = config%getFloats('tau0_twin',requiredSize=size(N_tw))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 16:23:36 +05:30
|
|
|
prm%c_1 = config%getFloat('twin_c',defaultVal=0.0_pReal)
|
|
|
|
prm%c_2 = config%getFloat('twin_b',defaultVal=1.0_pReal)
|
|
|
|
prm%c_3 = config%getFloat('twin_e',defaultVal=0.0_pReal)
|
|
|
|
prm%c_4 = config%getFloat('twin_d',defaultVal=0.0_pReal)
|
|
|
|
prm%gdot0_twin = config%getFloat('gdot0_twin')
|
|
|
|
prm%n_twin = config%getFloat('n_twin')
|
|
|
|
prm%spr = config%getFloat('s_pr')
|
|
|
|
prm%h0_TwinTwin = config%getFloat('h0_twintwin')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! expand: family => system
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_twin_0 = math_expand(xi_twin_0,N_tw)
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
! sanity checks
|
|
|
|
if (prm%gdot0_twin <= 0.0_pReal) extmsg = trim(extmsg)//' gdot0_twin'
|
|
|
|
if (prm%n_twin <= 0.0_pReal) extmsg = trim(extmsg)//' n_twin'
|
2020-03-16 16:23:36 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
else twinActive
|
2020-03-16 16:23:36 +05:30
|
|
|
xi_twin_0 = emptyRealArray
|
|
|
|
allocate(prm%gamma_twin_char,source=emptyRealArray)
|
2020-02-01 02:07:18 +05:30
|
|
|
allocate(prm%interaction_TwinTwin(0,0))
|
|
|
|
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
|
2020-02-01 02:07:18 +05:30
|
|
|
prm%h0_TwinSlip = config%getFloat('h0_twinslip')
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%interaction_SlipTwin = lattice_interaction_SlipByTwin(N_sl,N_tw,&
|
2020-02-01 02:07:18 +05:30
|
|
|
config%getFloats('interaction_sliptwin'), &
|
|
|
|
config%getString('lattice_structure'))
|
2020-03-16 03:37:23 +05:30
|
|
|
prm%interaction_TwinSlip = lattice_interaction_TwinBySlip(N_tw,N_sl,&
|
2020-02-01 02:07:18 +05:30
|
|
|
config%getFloats('interaction_twinslip'), &
|
|
|
|
config%getString('lattice_structure'))
|
|
|
|
else slipAndTwinActive
|
2020-03-16 01:46:28 +05:30
|
|
|
allocate(prm%interaction_SlipTwin(prm%sum_N_sl,prm%sum_N_tw)) ! at least one dimension is 0
|
|
|
|
allocate(prm%interaction_TwinSlip(prm%sum_N_tw,prm%sum_N_sl)) ! at least one dimension is 0
|
2020-02-01 02:07:18 +05:30
|
|
|
prm%h0_TwinSlip = 0.0_pReal
|
|
|
|
endif slipAndTwinActive
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-10-14 12:56:32 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! output pararameters
|
2020-02-14 11:47:30 +05:30
|
|
|
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
|
|
|
|
|
2014-07-02 17:57:39 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! allocate state arrays
|
2020-02-01 02:07:18 +05:30
|
|
|
NipcMyPhase = count(material_phaseAt == p) * discretization_nIP
|
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-04-01 13:26:59 +05:30
|
|
|
call material_allocateState(plasticState(p),NipcMyPhase,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
|
2020-02-01 02:07:18 +05:30
|
|
|
stt%xi_slip => plasticState(p)%state (startIndex:endIndex,:)
|
2020-03-16 16:23:36 +05:30
|
|
|
stt%xi_slip = spread(xi_slip_0, 2, NipcMyPhase)
|
2020-02-01 02:07:18 +05:30
|
|
|
dot%xi_slip => plasticState(p)%dotState(startIndex:endIndex,:)
|
2020-03-16 04:57:52 +05:30
|
|
|
plasticState(p)%atol(startIndex:endIndex) = config%getFloat('atol_xi',defaultVal=1.0_pReal)
|
|
|
|
if(any(plasticState(p)%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
|
2020-02-01 02:07:18 +05:30
|
|
|
stt%xi_twin => plasticState(p)%state (startIndex:endIndex,:)
|
2020-03-16 16:23:36 +05:30
|
|
|
stt%xi_twin = spread(xi_twin_0, 2, NipcMyPhase)
|
2020-02-01 02:07:18 +05:30
|
|
|
dot%xi_twin => plasticState(p)%dotState(startIndex:endIndex,:)
|
2020-03-16 04:57:52 +05:30
|
|
|
plasticState(p)%atol(startIndex:endIndex) = config%getFloat('atol_xi',defaultVal=1.0_pReal)
|
|
|
|
if(any(plasticState(p)%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
|
2020-02-01 02:07:18 +05:30
|
|
|
stt%gamma_slip => plasticState(p)%state (startIndex:endIndex,:)
|
|
|
|
dot%gamma_slip => plasticState(p)%dotState(startIndex:endIndex,:)
|
2020-03-16 04:57:52 +05:30
|
|
|
plasticState(p)%atol(startIndex:endIndex) = config%getFloat('atol_gamma',defaultVal=1.0e-6_pReal)
|
|
|
|
if(any(plasticState(p)%atol(startIndex:endIndex) < 0.0_pReal)) extmsg = trim(extmsg)//' atol_gamma'
|
2020-02-01 02:07:18 +05:30
|
|
|
! global alias
|
2020-03-16 04:57:52 +05:30
|
|
|
plasticState(p)%slipRate => plasticState(p)%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
|
2020-02-01 02:07:18 +05:30
|
|
|
stt%gamma_twin => plasticState(p)%state (startIndex:endIndex,:)
|
|
|
|
dot%gamma_twin => plasticState(p)%dotState(startIndex:endIndex,:)
|
2020-03-16 03:37:23 +05:30
|
|
|
plasticState(p)%atol(startIndex:endIndex) = config%getFloat('atol_gamma',defaultVal=1.0e-6_pReal)
|
2020-03-16 04:57:52 +05:30
|
|
|
if(any(plasticState(p)%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
|
|
|
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
|
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
|
|
|
|
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'('//PLASTICITY_PHENOPOWERLAW_LABEL//')')
|
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
enddo
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2018-04-25 23:11:18 +05:30
|
|
|
end subroutine 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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-03 01:13:02 +05:30
|
|
|
pure module subroutine plastic_phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
|
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) :: &
|
|
|
|
instance, &
|
|
|
|
of
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
integer :: &
|
|
|
|
i,k,l,m,n
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_slip_pos,gdot_slip_neg, &
|
|
|
|
dgdot_dtauslip_pos,dgdot_dtauslip_neg
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%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
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
associate(prm => param(instance))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
call kinetics_slip(Mp,instance,of,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
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
call kinetics_twin(Mp,instance,of,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
|
|
|
|
2014-12-08 21:25:30 +05:30
|
|
|
end subroutine plastic_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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-03 01:13:02 +05:30
|
|
|
module subroutine plastic_phenopowerlaw_dotState(Mp,instance,of)
|
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) :: &
|
|
|
|
instance, &
|
|
|
|
of
|
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
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%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
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
associate(prm => param(instance), stt => state(instance), dot => dotState(instance))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
sumGamma = sum(stt%gamma_slip(:,of))
|
|
|
|
sumF = sum(stt%gamma_twin(:,of)/prm%gamma_twin_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-02-01 02:13:12 +05:30
|
|
|
c_SlipSlip = prm%h0_slipslip * (1.0_pReal + prm%c_1*sumF** prm%c_2)
|
|
|
|
c_TwinSlip = prm%h0_TwinSlip * sumGamma**prm%c_3
|
|
|
|
c_TwinTwin = prm%h0_TwinTwin * 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-02-01 02:07:18 +05:30
|
|
|
left_SlipSlip = 1.0_pReal + prm%H_int
|
|
|
|
xi_slip_sat_offset = prm%spr*sqrt(sumF)
|
|
|
|
right_SlipSlip = abs(1.0_pReal-stt%xi_slip(:,of) / (prm%xi_slip_sat+xi_slip_sat_offset)) **prm%a_slip &
|
|
|
|
* sign(1.0_pReal,1.0_pReal-stt%xi_slip(:,of) / (prm%xi_slip_sat+xi_slip_sat_offset))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2018-09-12 13:29:09 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! shear rates
|
2020-02-01 02:07:18 +05:30
|
|
|
call kinetics_slip(Mp,instance,of,gdot_slip_pos,gdot_slip_neg)
|
|
|
|
dot%gamma_slip(:,of) = abs(gdot_slip_pos+gdot_slip_neg)
|
|
|
|
call kinetics_twin(Mp,instance,of,dot%gamma_twin(:,of))
|
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
|
2020-02-01 02:07:18 +05:30
|
|
|
dot%xi_slip(:,of) = c_SlipSlip * left_SlipSlip * &
|
|
|
|
matmul(prm%interaction_SlipSlip,dot%gamma_slip(:,of)*right_SlipSlip) &
|
|
|
|
+ matmul(prm%interaction_SlipTwin,dot%gamma_twin(:,of))
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
dot%xi_twin(:,of) = c_TwinSlip * matmul(prm%interaction_TwinSlip,dot%gamma_slip(:,of)) &
|
|
|
|
+ c_TwinTwin * matmul(prm%interaction_TwinTwin,dot%gamma_twin(:,of))
|
|
|
|
end associate
|
2018-08-26 00:43:32 +05:30
|
|
|
|
2014-12-08 21:25:30 +05:30
|
|
|
end subroutine plastic_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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-03 01:13:02 +05:30
|
|
|
module subroutine plastic_phenopowerlaw_results(instance,group)
|
2019-04-04 11:22:36 +05:30
|
|
|
|
|
|
|
integer, intent(in) :: instance
|
|
|
|
character(len=*), intent(in) :: group
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2019-04-04 11:22:36 +05:30
|
|
|
integer :: o
|
|
|
|
|
|
|
|
associate(prm => param(instance), stt => state(instance))
|
2020-02-14 11:47:30 +05:30
|
|
|
outputsLoop: do o = 1,size(prm%output)
|
|
|
|
select case(trim(prm%output(o)))
|
|
|
|
|
|
|
|
case('resistance_slip')
|
2020-03-16 01:46:28 +05:30
|
|
|
if(prm%sum_N_sl>0) call results_writeDataset(group,stt%xi_slip, 'xi_sl', &
|
|
|
|
'resistance against plastic slip','Pa')
|
2020-02-14 11:47:30 +05:30
|
|
|
case('accumulatedshear_slip')
|
2020-03-16 01:46:28 +05:30
|
|
|
if(prm%sum_N_sl>0) call results_writeDataset(group,stt%gamma_slip,'gamma_sl', &
|
|
|
|
'plastic shear','1')
|
2020-02-14 11:47:30 +05:30
|
|
|
|
|
|
|
case('resistance_twin')
|
2020-03-16 01:46:28 +05:30
|
|
|
if(prm%sum_N_tw>0) call results_writeDataset(group,stt%xi_twin, 'xi_tw', &
|
|
|
|
'resistance against twinning','Pa')
|
2020-02-14 11:47:30 +05:30
|
|
|
case('accumulatedshear_twin')
|
2020-03-16 01:46:28 +05:30
|
|
|
if(prm%sum_N_tw>0) call results_writeDataset(group,stt%gamma_twin,'gamma_tw', &
|
|
|
|
'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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-12-21 19:26:32 +05:30
|
|
|
pure subroutine kinetics_slip(Mp,instance,of, &
|
|
|
|
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) :: &
|
|
|
|
instance, &
|
|
|
|
of
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), intent(out), dimension(param(instance)%sum_N_sl) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_slip_pos, &
|
|
|
|
gdot_slip_neg
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), intent(out), optional, dimension(param(instance)%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
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%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
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
associate(prm => param(instance), stt => state(instance))
|
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-03-16 16:23:36 +05:30
|
|
|
gdot_slip_pos = prm%gdot0_slip * merge(0.5_pReal,1.0_pReal, prm%nonSchmidActive) & ! 1/2 if non-Schmid active
|
2020-02-01 02:07:18 +05:30
|
|
|
* sign(abs(tau_slip_pos/stt%xi_slip(:,of))**prm%n_slip, tau_slip_pos)
|
|
|
|
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))
|
|
|
|
gdot_slip_neg = prm%gdot0_slip * 0.5_pReal & ! only used if non-Schmid active, always 1/2
|
|
|
|
* sign(abs(tau_slip_neg/stt%xi_slip(:,of))**prm%n_slip, tau_slip_neg)
|
|
|
|
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))
|
|
|
|
dgdot_dtau_slip_pos = gdot_slip_pos*prm%n_slip/tau_slip_pos
|
|
|
|
else where
|
|
|
|
dgdot_dtau_slip_pos = 0.0_pReal
|
|
|
|
end where
|
|
|
|
endif
|
|
|
|
if (present(dgdot_dtau_slip_neg)) then
|
|
|
|
where(dNeq0(gdot_slip_neg))
|
|
|
|
dgdot_dtau_slip_neg = gdot_slip_neg*prm%n_slip/tau_slip_neg
|
|
|
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-12-21 19:26:32 +05:30
|
|
|
pure subroutine kinetics_twin(Mp,instance,of,&
|
|
|
|
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) :: &
|
|
|
|
instance, &
|
|
|
|
of
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%sum_N_tw), intent(out) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
gdot_twin
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%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
|
|
|
|
2020-03-16 01:46:28 +05:30
|
|
|
real(pReal), dimension(param(instance)%sum_N_tw) :: &
|
2020-02-01 02:07:18 +05:30
|
|
|
tau_twin
|
|
|
|
integer :: i
|
2020-02-14 11:47:30 +05:30
|
|
|
|
2020-02-01 02:07:18 +05:30
|
|
|
associate(prm => param(instance), stt => state(instance))
|
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)
|
|
|
|
gdot_twin = (1.0_pReal-sum(stt%gamma_twin(:,of)/prm%gamma_twin_char)) & ! only twin in untwinned volume fraction
|
|
|
|
* prm%gdot0_twin*(abs(tau_twin)/stt%xi_twin(:,of))**prm%n_twin
|
|
|
|
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))
|
|
|
|
dgdot_dtau_twin = gdot_twin*prm%n_twin/tau_twin
|
|
|
|
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
|
|
|
|
2019-12-03 01:13:02 +05:30
|
|
|
end submodule plastic_phenopowerlaw
|