DAMASK_EICMD/src/phase_mechanical_plastic_is...

317 lines
13 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
2018-12-30 22:41:03 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
2018-12-30 16:03:27 +05:30
!> @brief material subroutine for isotropic plasticity
!> @details Isotropic Plasticity which resembles the phenopowerlaw plasticity without
!! resolving the stress on the slip systems. Will give the response of phenopowerlaw for an
!! untextured polycrystal
!--------------------------------------------------------------------------------------------------
submodule(phase:plastic) isotropic
type :: tParameters
real(pREAL) :: &
2019-03-28 03:39:45 +05:30
M, & !< Taylor factor
2019-03-28 11:09:07 +05:30
dot_gamma_0, & !< reference strain rate
2019-03-28 03:12:02 +05:30
n, & !< stress exponent
h_0, &
h, & !< hardening pre-factor
2019-03-28 11:09:07 +05:30
h_ln, &
2019-03-28 03:39:45 +05:30
xi_inf, & !< maximum critical stress
2019-03-28 03:12:02 +05:30
a, &
2019-03-28 03:39:45 +05:30
c_1, &
c_4, &
c_3, &
2020-03-15 00:33:57 +05:30
c_2
2019-03-28 03:12:02 +05:30
logical :: &
dilatation
2023-06-04 10:47:38 +05:30
character(len=pSTRLEN), allocatable, dimension(:) :: &
output
2019-03-28 03:12:02 +05:30
end type tParameters
type :: tIsotropicState
real(pREAL), pointer, dimension(:) :: &
xi
2019-03-28 03:12:02 +05:30
end type tIsotropicState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
2021-05-24 20:50:33 +05:30
type(tParameters), allocatable, dimension(:) :: param
type(tIsotropicState), allocatable, dimension(:) :: state
2018-12-30 22:41:03 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief Perform module initialization.
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
module function plastic_isotropic_init() result(myPlasticity)
2020-08-15 19:32:10 +05:30
logical, dimension(:), allocatable :: myPlasticity
2019-03-28 03:12:02 +05:30
integer :: &
2021-02-16 20:36:09 +05:30
ph, &
2021-03-05 01:46:36 +05:30
Nmembers, &
2020-07-02 00:52:05 +05:30
sizeState, sizeDotState
real(pREAL) :: &
xi_0 !< initial critical stress
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('isotropic')
2022-12-07 22:59:03 +05:30
if (count(myPlasticity) == 0) return
2020-07-02 00:52:05 +05:30
print'(/,1x,a)', '<<<+- phase:mechanical:plastic:isotropic init -+>>>'
print'(/,a,i0)', ' # phases: ',count(myPlasticity); flush(IO_STDOUT)
2021-02-13 17:37:35 +05:30
2022-05-28 00:23:16 +05:30
print'(/,1x,a)', 'T. Maiti and P. Eisenlohr, Scripta Materialia 145:3740, 2018'
print'( 1x,a)', 'https://doi.org/10.1016/j.scriptamat.2017.09.047'
2023-01-18 23:20:01 +05:30
phases => config_material%get_dict('phase')
allocate(param(phases%length))
allocate(state(phases%length))
extmsg = ''
2021-02-16 20:36:09 +05:30
do ph = 1, phases%length
2022-12-07 22:59:03 +05:30
if (.not. myPlasticity(ph)) cycle
2021-02-16 20:36:09 +05:30
associate(prm => param(ph), &
stt => state(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,1x,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)
2020-08-15 19:32:10 +05:30
#else
2023-06-04 10:47:38 +05:30
prm%output = pl%get_as1dStr('output',defaultVal=emptyStrArray)
2020-08-15 19:32:10 +05:30
#endif
2020-03-15 00:33:57 +05:30
xi_0 = pl%get_asReal('xi_0')
prm%xi_inf = pl%get_asReal('xi_inf')
prm%dot_gamma_0 = pl%get_asReal('dot_gamma_0')
prm%n = pl%get_asReal('n')
prm%h_0 = pl%get_asReal('h_0')
prm%h = pl%get_asReal('h', defaultVal=3.0_pREAL) ! match for fcc random polycrystal
prm%M = pl%get_asReal('M')
prm%h_ln = pl%get_asReal('h_ln', defaultVal=0.0_pREAL)
prm%c_1 = pl%get_asReal('c_1', defaultVal=0.0_pREAL)
prm%c_4 = pl%get_asReal('c_4', defaultVal=0.0_pREAL)
prm%c_3 = pl%get_asReal('c_3', defaultVal=0.0_pREAL)
prm%c_2 = pl%get_asReal('c_2', defaultVal=0.0_pREAL)
prm%a = pl%get_asReal('a')
2020-08-15 19:32:10 +05:30
prm%dilatation = pl%get_asBool('dilatation',defaultVal = .false.)
2018-12-30 22:41:03 +05:30
!--------------------------------------------------------------------------------------------------
! sanity checks
if (xi_0 < 0.0_pREAL) extmsg = trim(extmsg)//' xi_0'
if (prm%dot_gamma_0 <= 0.0_pREAL) extmsg = trim(extmsg)//' dot_gamma_0'
if (prm%n <= 0.0_pREAL) extmsg = trim(extmsg)//' n'
if (prm%a <= 0.0_pREAL) extmsg = trim(extmsg)//' a'
if (prm%M <= 0.0_pREAL) extmsg = trim(extmsg)//' M'
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2023-01-23 13:01:59 +05:30
Nmembers = count(material_ID_phase == ph)
sizeDotState = size(['xi'])
2019-03-28 03:12:02 +05:30
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
stt%xi => plasticState(ph)%state(1,:)
stt%xi = xi_0
plasticState(ph)%atol(1) = pl%get_asReal('atol_xi',defaultVal=1.0_pREAL)
if (plasticState(ph)%atol(1) < 0.0_pREAL) extmsg = trim(extmsg)//' atol_xi'
2021-02-16 20:36:09 +05:30
2019-03-28 03:12:02 +05:30
end associate
2020-03-15 00:33:57 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg))
2020-03-15 00:33:57 +05:30
end do
2020-08-15 19:32:10 +05:30
end function plastic_isotropic_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
2021-04-25 11:36:52 +05:30
module subroutine isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2019-03-28 03:12:02 +05:30
Lp !< plastic velocity gradient
real(pREAL), dimension(3,3,3,3), intent(out) :: &
2019-03-28 03:12:02 +05:30
dLp_dMp !< derivative of Lp with respect to the Mandel stress
real(pREAL), dimension(3,3), intent(in) :: &
2019-03-28 03:12:02 +05:30
Mp !< Mandel stress
integer, intent(in) :: &
ph, &
2021-04-25 11:36:52 +05:30
en
2020-07-03 20:15:11 +05:30
real(pREAL), dimension(3,3) :: &
2019-03-28 03:12:02 +05:30
Mp_dev !< deviatoric part of the Mandel stress
real(pREAL) :: &
2019-03-28 11:09:07 +05:30
dot_gamma, & !< strainrate
2019-03-28 03:12:02 +05:30
norm_Mp_dev, & !< norm of the deviatoric part of the Mandel stress
squarenorm_Mp_dev !< square of the norm of the deviatoric part of the Mandel stress
integer :: &
k, l, m, n
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-23 03:39:51 +05:30
Mp_dev = math_deviatoric33(Mp)
squarenorm_Mp_dev = math_tensordot(Mp_dev,Mp_dev)
norm_Mp_dev = sqrt(squarenorm_Mp_dev)
if (norm_Mp_dev > 0.0_pREAL) then
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pREAL) * norm_Mp_dev/(prm%M*stt%xi(en)))**prm%n
2021-07-23 03:39:51 +05:30
Lp = dot_gamma * Mp_dev/norm_Mp_dev
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLp_dMp(k,l,m,n) = (prm%n-1.0_pREAL) * Mp_dev(k,l)*Mp_dev(m,n) / squarenorm_Mp_dev
2021-07-23 03:39:51 +05:30
forall (k=1:3,l=1:3) &
dLp_dMp(k,l,k,l) = dLp_dMp(k,l,k,l) + 1.0_pREAL
2021-07-23 03:39:51 +05:30
forall (k=1:3,m=1:3) &
dLp_dMp(k,k,m,m) = dLp_dMp(k,k,m,m) - 1.0_pREAL/3.0_pREAL
2021-07-23 03:39:51 +05:30
dLp_dMp = dot_gamma * dLp_dMp / norm_Mp_dev
else
Lp = 0.0_pREAL
dLp_dMp = 0.0_pREAL
2021-07-23 03:39:51 +05:30
end if
2019-03-28 03:12:02 +05:30
end associate
2019-01-06 04:11:13 +05:30
2021-01-26 12:03:04 +05:30
end subroutine isotropic_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief Calculate inelastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
2021-04-13 00:51:15 +05:30
module subroutine plastic_isotropic_LiAndItsTangent(Li,dLi_dMi,Mi,ph,en)
real(pREAL), dimension(3,3), intent(out) :: &
2019-03-28 03:12:02 +05:30
Li !< inleastic velocity gradient
real(pREAL), dimension(3,3,3,3), intent(out) :: &
dLi_dMi !< derivative of Li with respect to Mandel stress
real(pREAL), dimension(3,3), intent(in) :: &
Mi !< Mandel stress
2020-03-16 18:03:14 +05:30
integer, intent(in) :: &
ph, &
2021-04-13 00:51:15 +05:30
en
2020-07-03 20:15:11 +05:30
real(pREAL) :: &
tr !< trace of spherical part of Mandel stress (= 3 x pressure)
2019-03-28 03:12:02 +05:30
integer :: &
k, l, m, n
2021-07-23 03:39:51 +05:30
associate(prm => param(ph), stt => state(ph))
2021-07-23 03:39:51 +05:30
tr=math_trace33(math_spherical33(Mi))
if (prm%dilatation .and. abs(tr) > 0.0_pREAL) then ! no stress or J2 plasticity --> Li and its derivative are zero
2021-07-23 03:39:51 +05:30
Li = math_I3 &
* prm%dot_gamma_0 * (3.0_pREAL*prm%M*stt%xi(en))**(-prm%n) &
* tr * abs(tr)**(prm%n-1.0_pREAL)
2021-07-23 03:39:51 +05:30
forall (k=1:3,l=1:3,m=1:3,n=1:3) dLi_dMi(k,l,m,n) = prm%n / tr * Li(k,l) * math_I3(m,n)
else
Li = 0.0_pREAL
dLi_dMi = 0.0_pREAL
end if
2019-03-28 03:12:02 +05:30
end associate
2019-01-06 04:11:13 +05:30
end subroutine plastic_isotropic_LiAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief Calculate the rate of change of microstructure.
!--------------------------------------------------------------------------------------------------
module function isotropic_dotState(Mp,ph,en) result(dotState)
real(pREAL), dimension(3,3), intent(in) :: &
2019-03-28 03:12:02 +05:30
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) :: &
2019-03-28 11:09:07 +05:30
dot_gamma, & !< strainrate
2019-03-28 03:39:45 +05:30
xi_inf_star, & !< saturation xi
2019-03-28 03:12:02 +05:30
norm_Mp !< norm of the (deviatoric) Mandel stress
associate(prm => param(ph), stt => state(ph), dot_xi => dotState(1))
norm_Mp = merge(sqrt(math_tensordot(Mp,Mp)), &
sqrt(math_tensordot(math_deviatoric33(Mp),math_deviatoric33(Mp))), &
prm%dilatation)
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pREAL) * norm_Mp /(prm%M*stt%xi(en))) **prm%n
if (dot_gamma > 1e-12_pREAL) then
2019-03-28 03:39:45 +05:30
if (dEq0(prm%c_1)) then
xi_inf_star = prm%xi_inf
2019-03-28 03:12:02 +05:30
else
2019-03-28 03:39:45 +05:30
xi_inf_star = prm%xi_inf &
+ asinh( (dot_gamma / prm%c_1)**(1.0_pREAL / prm%c_2))**(1.0_pREAL / prm%c_3) &
/ prm%c_4 * (dot_gamma / prm%dot_gamma_0)**(1.0_pREAL / prm%n)
end if
dot_xi = dot_gamma &
* ( prm%h_0 + prm%h_ln * log(dot_gamma) ) &
* sign(abs(1.0_pREAL - stt%xi(en)/xi_inf_star)**prm%a *prm%h, 1.0_pREAL-stt%xi(en)/xi_inf_star)
2019-03-28 03:12:02 +05:30
else
dot_xi = 0.0_pREAL
end if
2019-03-28 03:12:02 +05:30
end associate
end function isotropic_dotState
2018-12-30 19:07:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
module subroutine plastic_isotropic_result(ph,group)
integer, intent(in) :: ph
2019-04-04 11:22:36 +05:30
character(len=*), intent(in) :: group
integer :: o
associate(prm => param(ph), stt => state(ph))
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
case ('xi')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt%xi,group,trim(prm%output(o)), &
'resistance against plastic flow','Pa')
end select
end do outputsLoop
end associate
2023-01-19 22:07:45 +05:30
end subroutine plastic_isotropic_result
2021-01-26 05:41:32 +05:30
end submodule isotropic