DAMASK_EICMD/src/phase_mechanics_plastic_iso...

335 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
!--------------------------------------------------------------------------------------------------
2021-01-26 05:41:32 +05:30
submodule(constitutive:plastic) isotropic
type :: tParameters
2019-03-28 03:12:02 +05:30
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, &
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
integer :: &
of_debug = 0
logical :: &
dilatation
character(len=pStringLen), allocatable, dimension(:) :: &
output
2019-03-28 03:12:02 +05:30
end type tParameters
type :: tIsotropicState
2019-03-28 03:12:02 +05:30
real(pReal), pointer, dimension(:) :: &
2019-03-28 03:39:45 +05:30
xi, &
2019-03-28 11:09:07 +05:30
gamma
2019-03-28 03:12:02 +05:30
end type tIsotropicState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:) :: param
type(tIsotropicState), allocatable, dimension(:) :: &
2018-12-30 16:03:27 +05:30
dotState, &
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 :: &
Ninstances, &
p, &
2020-08-15 19:32:10 +05:30
i, &
Nconstituents, &
2020-07-02 00:52:05 +05:30
sizeState, sizeDotState
real(pReal) :: &
xi_0 !< initial critical stress
2019-03-28 03:12:02 +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
print'(/,a)', ' <<<+- plastic_isotropic init -+>>>'
2020-08-15 19:32:10 +05:30
myPlasticity = plastic_active('isotropic')
Ninstances = count(myPlasticity)
print'(a,i2)', ' # instances: ',Ninstances; flush(IO_STDOUT)
if(Ninstances == 0) return
2020-07-02 00:52:05 +05:30
print*, 'Maiti and Eisenlohr, Scripta Materialia 145:3740, 2018'
print*, 'https://doi.org/10.1016/j.scriptamat.2017.09.047'
allocate(param(Ninstances))
allocate(state(Ninstances))
allocate(dotState(Ninstances))
phases => config_material%get('phase')
2020-08-15 19:32:10 +05:30
i = 0
do p = 1, phases%length
phase => phases%get(p)
2020-11-18 01:54:40 +05:30
mech => phase%get('mechanics')
2020-08-15 19:32:10 +05:30
if(.not. myPlasticity(p)) cycle
i = i + 1
associate(prm => param(i), &
dot => dotState(i), &
stt => state(i))
2020-11-03 03:16:46 +05:30
pl => mech%get('plasticity')
2020-08-15 19:32:10 +05:30
#if defined (__GFORTRAN__)
2020-08-15 19:32:10 +05:30
prm%output = output_asStrings(pl)
#else
prm%output = pl%get_asStrings('output',defaultVal=emptyStringArray)
#endif
2020-03-15 00:33:57 +05:30
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
if (p==material_phaseAt(debugConstitutive%grain,debugConstitutive%element)) &
prm%of_debug = material_phasememberAt(debugConstitutive%grain,debugConstitutive%ip,debugConstitutive%element)
2018-12-30 20:39:51 +05:30
#endif
2020-08-15 19:32:10 +05:30
xi_0 = pl%get_asFloat('xi_0')
prm%xi_inf = pl%get_asFloat('xi_inf')
prm%dot_gamma_0 = pl%get_asFloat('dot_gamma_0')
prm%n = pl%get_asFloat('n')
prm%h_0 = pl%get_asFloat('h_0')
2020-08-15 19:32:10 +05:30
prm%M = pl%get_asFloat('M')
prm%h_ln = pl%get_asFloat('h_ln', defaultVal=0.0_pReal)
prm%c_1 = pl%get_asFloat('c_1', defaultVal=0.0_pReal)
prm%c_4 = pl%get_asFloat('c_4', defaultVal=0.0_pReal)
prm%c_3 = pl%get_asFloat('c_3', defaultVal=0.0_pReal)
prm%c_2 = pl%get_asFloat('c_2', defaultVal=0.0_pReal)
prm%a = pl%get_asFloat('a')
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
2021-01-26 12:24:24 +05:30
Nconstituents = count(material_phaseAt2 == p)
2020-08-15 19:32:10 +05:30
sizeDotState = size(['xi ','gamma'])
2019-03-28 03:12:02 +05:30
sizeState = sizeDotState
call constitutive_allocateState(plasticState(p),Nconstituents,sizeState,sizeDotState,0)
!--------------------------------------------------------------------------------------------------
2020-03-16 03:37:23 +05:30
! state aliases and initialization
2019-03-28 03:39:45 +05:30
stt%xi => plasticState(p)%state (1,:)
stt%xi = xi_0
2019-03-28 03:39:45 +05:30
dot%xi => plasticState(p)%dotState(1,:)
2020-08-15 19:32:10 +05:30
plasticState(p)%atol(1) = pl%get_asFloat('atol_xi',defaultVal=1.0_pReal)
if (plasticState(p)%atol(1) < 0.0_pReal) extmsg = trim(extmsg)//' atol_xi'
2019-03-28 11:09:07 +05:30
stt%gamma => plasticState(p)%state (2,:)
dot%gamma => plasticState(p)%dotState(2,:)
2020-08-15 19:32:10 +05:30
plasticState(p)%atol(2) = pl%get_asFloat('atol_gamma',defaultVal=1.0e-6_pReal)
if (plasticState(p)%atol(2) < 0.0_pReal) extmsg = trim(extmsg)//' atol_gamma'
2019-03-28 03:12:02 +05:30
! global alias
plasticState(p)%slipRate => plasticState(p)%dotState(2:2,:)
2019-03-28 03:12:02 +05:30
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
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
2020-08-15 19:32:10 +05:30
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(isotropic)')
2020-03-15 00:33:57 +05:30
2019-03-28 03:12:02 +05:30
enddo
2020-08-15 19:32:10 +05:30
end function plastic_isotropic_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
2021-01-26 13:09:17 +05:30
module subroutine isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,me)
2019-03-28 03:12:02 +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
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
instance, &
2021-01-26 13:09:17 +05:30
me
2020-07-03 20:15:11 +05:30
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3) :: &
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
2019-03-28 03:12:02 +05:30
associate(prm => param(instance), stt => state(instance))
2019-03-28 03:12:02 +05:30
Mp_dev = math_deviatoric33(Mp)
2020-03-16 23:13:04 +05:30
squarenorm_Mp_dev = math_tensordot(Mp_dev,Mp_dev)
2019-03-28 03:12:02 +05:30
norm_Mp_dev = sqrt(squarenorm_Mp_dev)
2019-03-28 03:12:02 +05:30
if (norm_Mp_dev > 0.0_pReal) then
2021-01-26 13:09:17 +05:30
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pReal) * norm_Mp_dev/(prm%M*stt%xi(me))) **prm%n
Lp = dot_gamma/prm%M * Mp_dev/norm_Mp_dev
2019-03-28 03:12:02 +05:30
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
forall (k=1:3,l=1:3) &
dLp_dMp(k,l,k,l) = dLp_dMp(k,l,k,l) + 1.0_pReal
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
2019-03-28 11:09:07 +05:30
dLp_dMp = dot_gamma / prm%M * dLp_dMp / norm_Mp_dev
2019-03-28 03:12:02 +05:30
else
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
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-01-26 13:09:17 +05:30
module subroutine plastic_isotropic_LiAndItsTangent(Li,dLi_dMi,Mi,instance,me)
2020-03-16 18:03:14 +05:30
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
2020-03-16 18:03:14 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mi !< Mandel stress
2020-03-16 18:03:14 +05:30
integer, intent(in) :: &
2019-03-28 03:12:02 +05:30
instance, &
2021-01-26 13:09:17 +05:30
me
2020-07-03 20:15:11 +05:30
2019-03-28 03:12:02 +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
2019-03-28 03:12:02 +05:30
associate(prm => param(instance), stt => state(instance))
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
Li = math_I3 &
2021-01-26 13:09:17 +05:30
* prm%dot_gamma_0/prm%M * (3.0_pReal*prm%M*stt%xi(me))**(-prm%n) &
* tr * abs(tr)**(prm%n-1.0_pReal)
2021-01-26 13:09:17 +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)
2019-03-28 03:12:02 +05:30
else
Li = 0.0_pReal
dLi_dMi = 0.0_pReal
2019-03-28 03:12:02 +05:30
endif
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.
!--------------------------------------------------------------------------------------------------
2021-01-26 16:47:00 +05:30
module subroutine isotropic_dotState(Mp,instance,me)
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
instance, &
2021-01-26 13:09:17 +05:30
me
2019-03-28 03:12:02 +05:30
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
2019-03-28 03:12:02 +05:30
associate(prm => param(instance), stt => state(instance), dot => dotState(instance))
2019-03-28 03:12:02 +05:30
if (prm%dilatation) then
2020-03-16 23:13:04 +05:30
norm_Mp = sqrt(math_tensordot(Mp,Mp))
2019-03-28 03:12:02 +05:30
else
2020-03-16 23:13:04 +05:30
norm_Mp = sqrt(math_tensordot(math_deviatoric33(Mp),math_deviatoric33(Mp)))
2019-03-28 03:12:02 +05:30
endif
2021-01-26 13:09:17 +05:30
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pReal) * norm_Mp /(prm%M*stt%xi(me))) **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 &
2019-04-04 11:22:36 +05:30
+ 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)
2019-03-28 03:12:02 +05:30
endif
2021-01-26 13:09:17 +05:30
dot%xi(me) = dot_gamma &
* ( prm%h_0 + prm%h_ln * log(dot_gamma) ) &
2021-01-26 13:09:17 +05:30
* abs( 1.0_pReal - stt%xi(me)/xi_inf_star )**prm%a &
* sign(1.0_pReal, 1.0_pReal - stt%xi(me)/xi_inf_star)
2019-03-28 03:12:02 +05:30
else
2021-01-26 13:09:17 +05:30
dot%xi(me) = 0.0_pReal
2019-03-28 03:12:02 +05:30
endif
2021-01-26 13:09:17 +05:30
dot%gamma(me) = dot_gamma ! ToDo: not really used
2019-03-28 03:12:02 +05:30
end associate
2021-01-26 16:47:00 +05:30
end subroutine isotropic_dotState
2018-12-30 19:07:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_isotropic_results(instance,group)
integer, intent(in) :: instance
2019-04-04 11:22:36 +05:30
character(len=*), intent(in) :: group
integer :: o
associate(prm => param(instance), stt => state(instance))
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
case ('xi')
2020-08-15 19:32:10 +05:30
call results_writeDataset(group,stt%xi,trim(prm%output(o)), &
'resistance against plastic flow','Pa')
end select
enddo outputsLoop
end associate
end subroutine plastic_isotropic_results
2021-01-26 05:41:32 +05:30
end submodule isotropic