DAMASK_EICMD/src/constitutive_plastic_isotro...

335 lines
14 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(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
xi_0, & !< initial critical stress
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
h0, &
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 subroutine plastic_isotropic_init
2019-03-28 03:12:02 +05:30
integer :: &
Ninstance, &
p, &
2019-03-28 03:12:02 +05:30
NipcMyPhase, &
sizeState, sizeDotState
2019-03-28 03:12:02 +05:30
character(len=pStringLen) :: &
extmsg = ''
2020-03-15 00:33:57 +05:30
write(6,'(/,a)') ' <<<+- plastic_'//PLASTICITY_ISOTROPIC_LABEL//' init -+>>>'; flush(6)
write(6,'(/,a)') ' Maiti and Eisenlohr, Scripta Materialia 145:3740, 2018'
write(6,'(a)') ' https://doi.org/10.1016/j.scriptamat.2017.09.047'
2019-03-28 03:12:02 +05:30
Ninstance = count(phase_plasticity == PLASTICITY_ISOTROPIC_ID)
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-03-28 03:12:02 +05:30
allocate(param(Ninstance))
allocate(state(Ninstance))
allocate(dotState(Ninstance))
2019-03-28 03:12:02 +05:30
do p = 1, size(phase_plasticity)
if (phase_plasticity(p) /= PLASTICITY_ISOTROPIC_ID) cycle
associate(prm => param(phase_plasticityInstance(p)), &
dot => dotState(phase_plasticityInstance(p)), &
stt => state(phase_plasticityInstance(p)), &
config => config_phase(p))
2020-03-15 00:33:57 +05:30
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
if (p==material_phaseAt(debug_g,debug_e)) &
2019-06-14 12:47:05 +05:30
prm%of_debug = material_phasememberAt(debug_g,debug_i,debug_e)
2018-12-30 20:39:51 +05:30
#endif
2019-03-28 03:39:45 +05:30
prm%xi_0 = config%getFloat('tau0')
prm%xi_inf = config%getFloat('tausat')
2019-03-28 11:09:07 +05:30
prm%dot_gamma_0 = config%getFloat('gdot0')
2019-03-28 03:12:02 +05:30
prm%n = config%getFloat('n')
prm%h0 = config%getFloat('h0')
2019-03-28 03:39:45 +05:30
prm%M = config%getFloat('m')
2019-03-28 11:09:07 +05:30
prm%h_ln = config%getFloat('h0_slopelnrate', defaultVal=0.0_pReal)
2019-03-28 03:39:45 +05:30
prm%c_1 = config%getFloat('tausat_sinhfita',defaultVal=0.0_pReal)
prm%c_4 = config%getFloat('tausat_sinhfitb',defaultVal=0.0_pReal)
prm%c_3 = config%getFloat('tausat_sinhfitc',defaultVal=0.0_pReal)
prm%c_2 = config%getFloat('tausat_sinhfitd',defaultVal=0.0_pReal)
2019-03-28 03:12:02 +05:30
prm%a = config%getFloat('a')
2019-03-28 03:12:02 +05:30
prm%dilatation = config%keyExists('/dilatation/')
2018-12-30 22:41:03 +05:30
!--------------------------------------------------------------------------------------------------
! sanity checks
if (prm%xi_0 < 0.0_pReal) extmsg = trim(extmsg)//' xi_0'
2019-03-28 11:09:07 +05:30
if (prm%dot_gamma_0 <= 0.0_pReal) extmsg = trim(extmsg)//' dot_gamma_0'
2019-03-28 03:12:02 +05:30
if (prm%n <= 0.0_pReal) extmsg = trim(extmsg)//' n'
if (prm%a <= 0.0_pReal) extmsg = trim(extmsg)//' a'
2019-03-28 03:39:45 +05:30
if (prm%M <= 0.0_pReal) extmsg = trim(extmsg)//' m'
!--------------------------------------------------------------------------------------------------
! allocate state arrays
NipcMyPhase = count(material_phaseAt == p) * discretization_nIP
2019-03-28 03:39:45 +05:30
sizeDotState = size(['xi ','accumulated_shear'])
2019-03-28 03:12:02 +05:30
sizeState = sizeDotState
call material_allocatePlasticState(p,NipcMyPhase,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 = prm%xi_0
dot%xi => plasticState(p)%dotState(1,:)
plasticState(p)%atol(1) = config%getFloat('atol_flowstress',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,:)
plasticState(p)%atol(2) = config%getFloat('atol_shear',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
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'('//PLASTICITY_ISOTROPIC_LABEL//')')
2019-03-28 03:12:02 +05:30
enddo
end subroutine plastic_isotropic_init
!--------------------------------------------------------------------------------------------------
!> @brief Calculate plastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
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, &
of
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)
squarenorm_Mp_dev = math_mul33xx33(Mp_dev,Mp_dev)
norm_Mp_dev = sqrt(squarenorm_Mp_dev)
2019-03-28 03:12:02 +05:30
if (norm_Mp_dev > 0.0_pReal) then
2019-03-28 11:09:07 +05:30
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pReal) * norm_Mp_dev/(prm%M*stt%xi(of))) **prm%n
Lp = dot_gamma/prm%M * Mp_dev/norm_Mp_dev
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
2019-03-28 03:12:02 +05:30
if (iand(debug_level(debug_constitutive), debug_levelExtensive) /= 0 &
.and. (of == prm%of_debug .or. .not. iand(debug_level(debug_constitutive),debug_levelSelective) /= 0)) then
write(6,'(/,a,/,3(12x,3(f12.4,1x)/))') '<< CONST isotropic >> Tstar (dev) / MPa', &
transpose(Mp_dev)*1.0e-6_pReal
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> norm Tstar / MPa', norm_Mp_dev*1.0e-6_pReal
2019-03-28 11:09:07 +05:30
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> gdot', dot_gamma
2019-03-28 03:12:02 +05:30
end if
2018-12-30 20:39:51 +05:30
#endif
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
end subroutine plastic_isotropic_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief Calculate inelastic velocity gradient and its tangent.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_isotropic_LiAndItsTangent(Li,dLi_dMi,Mi,instance,of)
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3), intent(out) :: &
Li !< inleastic velocity gradient
real(pReal), dimension(3,3,3,3), intent(out) :: &
dLi_dMi !< derivative of Li with respect to Mandel stress
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mi !< Mandel stress
2019-03-28 03:12:02 +05:30
integer, intent(in) :: &
instance, &
of
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 &
* prm%dot_gamma_0/prm%M * (3.0_pReal*prm%M*stt%xi(of))**(-prm%n) &
* tr * abs(tr)**(prm%n-1.0_pReal)
#ifdef DEBUG
if (iand(debug_level(debug_constitutive), debug_levelExtensive) /= 0 &
.and. (of == prm%of_debug .or. .not. iand(debug_level(debug_constitutive),debug_levelSelective) /= 0)) then
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> pressure / MPa', tr/3.0_pReal*1.0e-6_pReal
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> gdot', prm%dot_gamma_0 * (3.0_pReal*prm%M*stt%xi(of))**(-prm%n) &
* tr * abs(tr)**(prm%n-1.0_pReal)
end if
#endif
2019-03-28 03:12:02 +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.
!--------------------------------------------------------------------------------------------------
module subroutine plastic_isotropic_dotState(Mp,instance,of)
2019-03-28 03:12:02 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer, intent(in) :: &
instance, &
of
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
norm_Mp = sqrt(math_mul33xx33(Mp,Mp))
else
norm_Mp = sqrt(math_mul33xx33(math_deviatoric33(Mp),math_deviatoric33(Mp)))
endif
2019-03-28 11:09:07 +05:30
dot_gamma = prm%dot_gamma_0 * (sqrt(1.5_pReal) * norm_Mp /(prm%M*stt%xi(of))) **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
dot%xi(of) = dot_gamma &
* ( prm%h0 + prm%h_ln * log(dot_gamma) ) &
* abs( 1.0_pReal - stt%xi(of)/xi_inf_star )**prm%a &
* sign(1.0_pReal, 1.0_pReal - stt%xi(of)/xi_inf_star)
2019-03-28 03:12:02 +05:30
else
2019-03-28 11:09:07 +05:30
dot%xi(of) = 0.0_pReal
2019-03-28 03:12:02 +05:30
endif
2019-03-28 11:09:07 +05:30
dot%gamma(of) = dot_gamma ! ToDo: not really used
2019-03-28 03:12:02 +05:30
end associate
end subroutine plastic_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 ('flowstress') ! ToDo: should be 'xi'
2019-04-04 11:22:36 +05:30
call results_writeDataset(group,stt%xi,'xi','resistance against plastic flow','Pa')
end select
enddo outputsLoop
end associate
end subroutine plastic_isotropic_results
end submodule plastic_isotropic