DAMASK_EICMD/src/plastic_isotropic.f90

513 lines
20 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
!--------------------------------------------------------------------------------------------------
module plastic_isotropic
use prec, only: &
2019-01-06 04:11:13 +05:30
pReal, &
pInt
2018-12-30 22:41:03 +05:30
implicit none
private
integer(pInt), dimension(:,:), allocatable, target, public :: &
2018-12-30 16:03:27 +05:30
plastic_isotropic_sizePostResult !< size of each post result output
character(len=64), dimension(:,:), allocatable, target, public :: &
2018-12-30 16:03:27 +05:30
plastic_isotropic_output !< name of each post result output
2018-12-30 22:41:03 +05:30
enum, bind(c)
2018-12-30 16:03:27 +05:30
enumerator :: &
undefined_ID, &
flowstress_ID, &
strainrate_ID
end enum
2018-12-30 16:03:27 +05:30
type, private :: tParameters
2018-12-30 20:39:51 +05:30
real(pReal) :: &
2018-12-30 16:03:27 +05:30
fTaylor, & !< Taylor factor
tau0, & !< initial critical stress
gdot0, & !< reference strain rate
n, & !< stress exponent
h0, &
h0_slopeLnRate, &
2018-12-30 16:03:27 +05:30
tausat, & !< maximum critical stress
a, &
tausat_SinhFitA, &
tausat_SinhFitB, &
tausat_SinhFitC, &
2018-12-30 16:03:27 +05:30
tausat_SinhFitD, &
aTolFlowstress, &
aTolShear
2018-12-30 20:39:51 +05:30
integer(pInt) :: &
of_debug = 0_pInt
2019-01-07 11:54:02 +05:30
integer(kind(undefined_ID)), allocatable, dimension(:) :: &
2018-12-30 16:03:27 +05:30
outputID
2018-12-30 20:39:51 +05:30
logical :: &
dilatation
end type tParameters
2018-12-30 16:03:27 +05:30
type, private :: tIsotropicState
2018-12-30 22:41:03 +05:30
real(pReal), pointer, dimension(:) :: &
flowstress, &
accumulatedShear
end type tIsotropicState
!--------------------------------------------------------------------------------------------------
! containers for parameters and state
type(tParameters), allocatable, dimension(:), private :: param
2018-12-30 16:03:27 +05:30
type(tIsotropicState), allocatable, dimension(:), private :: &
dotState, &
state
2018-12-30 22:41:03 +05:30
public :: &
plastic_isotropic_init, &
plastic_isotropic_LpAndItsTangent, &
plastic_isotropic_LiAndItsTangent, &
plastic_isotropic_dotState, &
plastic_isotropic_postResults, &
plastic_isotropic_results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2019-03-09 05:37:26 +05:30
subroutine plastic_isotropic_init
2018-12-30 22:41:03 +05:30
use prec, only: &
pStringLen
use debug, only: &
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
debug_e, &
debug_i, &
debug_g, &
debug_levelExtensive, &
#endif
debug_level, &
2019-01-06 04:25:10 +05:30
debug_constitutive, &
debug_levelBasic
2018-12-30 16:03:27 +05:30
use IO, only: &
2019-03-09 05:37:26 +05:30
IO_error
use material, only: &
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
phasememberAt, &
#endif
phase_plasticity, &
phase_plasticityInstance, &
phase_Noutput, &
2018-12-30 16:03:27 +05:30
material_allocatePlasticState, &
PLASTICITY_ISOTROPIC_label, &
PLASTICITY_ISOTROPIC_ID, &
material_phase, &
plasticState
use config, only: &
2018-06-27 00:24:54 +05:30
config_phase
2018-12-30 22:41:03 +05:30
use lattice
implicit none
integer(pInt) :: &
2018-12-30 22:41:03 +05:30
Ninstance, &
p, i, &
NipcMyPhase, &
sizeState, sizeDotState
2019-01-07 11:54:02 +05:30
2018-12-30 22:41:03 +05:30
character(len=65536), dimension(0), parameter :: emptyStringArray = [character(len=65536)::]
2019-01-06 04:11:13 +05:30
2018-12-30 16:03:27 +05:30
integer(kind(undefined_ID)) :: &
2018-12-30 22:41:03 +05:30
outputID
2018-12-30 16:03:27 +05:30
2018-12-30 22:41:03 +05:30
character(len=pStringLen) :: &
extmsg = ''
character(len=65536), dimension(:), allocatable :: &
outputs
2018-12-30 22:41:03 +05:30
write(6,'(/,a)') ' <<<+- plastic_'//PLASTICITY_ISOTROPIC_label//' init -+>>>'
2019-03-09 15:32:12 +05:30
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-09 05:37:26 +05:30
Ninstance = count(phase_plasticity == PLASTICITY_ISOTROPIC_ID)
2018-12-30 22:41:03 +05:30
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0_pInt) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2018-12-30 20:39:51 +05:30
2019-01-07 12:34:02 +05:30
allocate(plastic_isotropic_sizePostResult(maxval(phase_Noutput),Ninstance),source=0_pInt)
allocate(plastic_isotropic_output(maxval(phase_Noutput),Ninstance))
plastic_isotropic_output = ''
2018-12-30 22:41:03 +05:30
allocate(param(Ninstance))
allocate(state(Ninstance))
allocate(dotState(Ninstance))
2019-01-06 04:25:10 +05:30
do p = 1_pInt, size(phase_plasticity)
2018-12-30 16:03:27 +05:30
if (phase_plasticity(p) /= PLASTICITY_ISOTROPIC_ID) cycle
2018-12-30 22:41:03 +05:30
associate(prm => param(phase_plasticityInstance(p)), &
dot => dotState(phase_plasticityInstance(p)), &
stt => state(phase_plasticityInstance(p)), &
config => config_phase(p))
2019-01-07 11:54:02 +05:30
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
if (p==material_phase(debug_g,debug_i,debug_e)) then
2018-12-30 22:41:03 +05:30
prm%of_debug = phasememberAt(debug_g,debug_i,debug_e)
2018-12-30 20:39:51 +05:30
endif
#endif
2019-01-07 12:34:02 +05:30
prm%tau0 = config%getFloat('tau0')
prm%tausat = config%getFloat('tausat')
prm%gdot0 = config%getFloat('gdot0')
prm%n = config%getFloat('n')
prm%h0 = config%getFloat('h0')
prm%fTaylor = config%getFloat('m')
prm%h0_slopeLnRate = config%getFloat('h0_slopelnrate', defaultVal=0.0_pReal)
prm%tausat_SinhFitA = config%getFloat('tausat_sinhfita',defaultVal=0.0_pReal)
prm%tausat_SinhFitB = config%getFloat('tausat_sinhfitb',defaultVal=0.0_pReal)
prm%tausat_SinhFitC = config%getFloat('tausat_sinhfitc',defaultVal=0.0_pReal)
prm%tausat_SinhFitD = config%getFloat('tausat_sinhfitd',defaultVal=0.0_pReal)
prm%a = config%getFloat('a')
prm%aTolFlowStress = config%getFloat('atol_flowstress',defaultVal=1.0_pReal)
prm%aTolShear = config%getFloat('atol_shear', defaultVal=1.0e-6_pReal)
2019-01-07 11:54:02 +05:30
2018-12-30 22:41:03 +05:30
prm%dilatation = config%keyExists('/dilatation/')
!--------------------------------------------------------------------------------------------------
! sanity checks
extmsg = ''
2019-01-07 12:34:02 +05:30
if (prm%aTolShear <= 0.0_pReal) extmsg = trim(extmsg)//' aTolShear'
if (prm%tau0 < 0.0_pReal) extmsg = trim(extmsg)//' tau0'
if (prm%gdot0 <= 0.0_pReal) extmsg = trim(extmsg)//' gdot0'
if (prm%n <= 0.0_pReal) extmsg = trim(extmsg)//' n'
if (prm%tausat <= prm%tau0) extmsg = trim(extmsg)//' tausat'
if (prm%a <= 0.0_pReal) extmsg = trim(extmsg)//' a'
if (prm%fTaylor <= 0.0_pReal) extmsg = trim(extmsg)//' m'
if (prm%aTolFlowstress <= 0.0_pReal) extmsg = trim(extmsg)//' atol_flowstress'
if (prm%aTolShear <= 0.0_pReal) extmsg = trim(extmsg)//' atol_shear'
2019-01-07 11:54:02 +05:30
2018-12-30 22:41:03 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') &
call IO_error(211_pInt,ext_msg=trim(extmsg)//'('//PLASTICITY_ISOTROPIC_label//')')
!--------------------------------------------------------------------------------------------------
! output pararameters
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
2018-12-30 16:03:27 +05:30
allocate(prm%outputID(0))
do i=1_pInt, size(outputs)
outputID = undefined_ID
select case(outputs(i))
2019-01-06 04:11:13 +05:30
2018-12-30 16:03:27 +05:30
case ('flowstress')
outputID = flowstress_ID
case ('strainrate')
outputID = strainrate_ID
2019-01-06 04:11:13 +05:30
2018-12-30 16:03:27 +05:30
end select
2019-01-06 04:11:13 +05:30
2018-12-30 16:03:27 +05:30
if (outputID /= undefined_ID) then
2019-01-06 04:11:13 +05:30
plastic_isotropic_output(i,phase_plasticityInstance(p)) = outputs(i)
plastic_isotropic_sizePostResult(i,phase_plasticityInstance(p)) = 1_pInt
prm%outputID = [prm%outputID, outputID]
2018-12-30 16:03:27 +05:30
endif
2019-01-06 04:25:10 +05:30
enddo
!--------------------------------------------------------------------------------------------------
! allocate state arrays
2018-12-30 22:41:03 +05:30
NipcMyPhase = count(material_phase == p)
2019-01-06 04:11:13 +05:30
sizeDotState = size(['flowstress ','accumulated_shear'])
sizeState = sizeDotState
2018-12-30 16:03:27 +05:30
call material_allocatePlasticState(p,NipcMyPhase,sizeState,sizeDotState,0_pInt, &
1_pInt,0_pInt,0_pInt)
plasticState(p)%sizePostResults = sum(plastic_isotropic_sizePostResult(:,phase_plasticityInstance(p)))
!--------------------------------------------------------------------------------------------------
! locally defined state aliases and initialization of state0 and aTolState
2019-01-07 12:34:02 +05:30
stt%flowstress => plasticState(p)%state (1,:)
2018-12-30 22:41:03 +05:30
stt%flowstress = prm%tau0
2019-01-07 12:34:02 +05:30
dot%flowstress => plasticState(p)%dotState(1,:)
plasticState(p)%aTolState(1) = prm%aTolFlowstress
2019-01-07 12:34:02 +05:30
stt%accumulatedShear => plasticState(p)%state (2,:)
dot%accumulatedShear => plasticState(p)%dotState(2,:)
plasticState(p)%aTolState(2) = prm%aTolShear
2018-12-30 16:03:27 +05:30
! global alias
2019-01-07 12:34:02 +05:30
plasticState(p)%slipRate => plasticState(p)%dotState(2:2,:)
plasticState(p)%accumulatedSlip => plasticState(p)%state (2:2,:)
2019-01-07 11:54:02 +05:30
2018-12-30 22:41:03 +05:30
plasticState(p)%state0 = plasticState(p)%state ! ToDo: this could be done centrally
2019-01-06 04:11:13 +05:30
end associate
enddo
end subroutine plastic_isotropic_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
!--------------------------------------------------------------------------------------------------
2018-12-30 20:39:51 +05:30
subroutine plastic_isotropic_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
#ifdef DEBUG
use debug, only: &
debug_level, &
2018-12-30 20:39:51 +05:30
debug_constitutive,&
debug_levelExtensive, &
2018-12-30 20:39:51 +05:30
debug_levelSelective
#endif
use math, only: &
math_deviatoric33, &
math_mul33xx33
implicit none
2018-12-30 17:05:26 +05:30
real(pReal), dimension(3,3), intent(out) :: &
Lp !< plastic velocity gradient
2018-12-30 17:05:26 +05:30
real(pReal), dimension(3,3,3,3), intent(out) :: &
dLp_dMp !< derivative of Lp with respect to the Mandel stress
2018-12-30 17:05:26 +05:30
real(pReal), dimension(3,3), intent(in) :: &
2018-12-30 20:39:51 +05:30
Mp !< Mandel stress
integer(pInt), intent(in) :: &
2018-12-30 20:39:51 +05:30
instance, &
of
real(pReal), dimension(3,3) :: &
2018-12-30 17:05:26 +05:30
Mp_dev !< deviatoric part of the Mandel stress
real(pReal) :: &
gamma_dot, & !< strainrate
2018-12-30 20:39:51 +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(pInt) :: &
k, l, m, n
2018-12-30 20:39:51 +05:30
associate(prm => param(instance), stt => state(instance))
2019-01-07 11:54:02 +05:30
2018-12-30 17:05:26 +05:30
Mp_dev = math_deviatoric33(Mp)
squarenorm_Mp_dev = math_mul33xx33(Mp_dev,Mp_dev)
2019-01-07 11:54:02 +05:30
norm_Mp_dev = sqrt(squarenorm_Mp_dev)
2018-12-30 20:39:51 +05:30
if (norm_Mp_dev > 0.0_pReal) then
2018-12-30 22:41:03 +05:30
gamma_dot = prm%gdot0 * (sqrt(1.5_pReal) * norm_Mp_dev/(prm%fTaylor*stt%flowstress(of))) **prm%n
2019-01-07 11:54:02 +05:30
Lp = Mp_dev/norm_Mp_dev * gamma_dot/prm%fTaylor
2018-12-30 20:39:51 +05:30
#ifdef DEBUG
if (iand(debug_level(debug_constitutive), debug_levelExtensive) /= 0_pInt &
2018-12-30 20:39:51 +05:30
.and. (of == prm%of_debug .or. .not. iand(debug_level(debug_constitutive),debug_levelSelective) /= 0_pInt)) then
write(6,'(/,a,/,3(12x,3(f12.4,1x)/))') '<< CONST isotropic >> Tstar (dev) / MPa', &
2018-12-30 17:05:26 +05:30
transpose(Mp_dev)*1.0e-6_pReal
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> norm Tstar / MPa', norm_Mp_dev*1.0e-6_pReal
write(6,'(/,a,/,f12.5)') '<< CONST isotropic >> gdot', gamma_dot
end if
2018-12-30 20:39:51 +05:30
#endif
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) &
2018-12-30 17:05:26 +05:30
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_pInt:3_pInt,l=1_pInt:3_pInt) &
2018-12-30 17:05:26 +05:30
dLp_dMp(k,l,k,l) = dLp_dMp(k,l,k,l) + 1.0_pReal
forall (k=1_pInt:3_pInt,m=1_pInt:3_pInt) &
2018-12-30 17:05:26 +05:30
dLp_dMp(k,k,m,m) = dLp_dMp(k,k,m,m) - 1.0_pReal/3.0_pReal
dLp_dMp = gamma_dot / prm%fTaylor * dLp_dMp / norm_Mp_dev
2018-12-30 20:39:51 +05:30
else
Lp = 0.0_pReal
dLp_dMp = 0.0_pReal
end if
2019-01-07 11:54:02 +05:30
2018-12-30 16:03:27 +05:30
end associate
2019-01-06 04:11:13 +05:30
end subroutine plastic_isotropic_LpAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief calculates plastic velocity gradient and its tangent
2018-12-30 22:41:03 +05:30
! ToDo: Rename Tstar to Mi?
!--------------------------------------------------------------------------------------------------
subroutine plastic_isotropic_LiAndItsTangent(Li,dLi_dTstar,Tstar,instance,of)
use math, only: &
math_spherical33, &
math_mul33xx33
implicit none
real(pReal), dimension(3,3), intent(out) :: &
Li !< inleastic velocity gradient
real(pReal), dimension(3,3,3,3), intent(out) :: &
dLi_dTstar !< derivative of Li with respect to the Mandel stress
2019-01-07 11:54:02 +05:30
2018-12-30 17:05:26 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Tstar !< Mandel stress ToDo: Mi?
2018-12-30 19:07:31 +05:30
integer(pInt), intent(in) :: &
instance, &
of
real(pReal), dimension(3,3) :: &
Tstar_sph !< sphiatoric part of the Mandel stress
real(pReal) :: &
gamma_dot, & !< strainrate
norm_Tstar_sph, & !< euclidean norm of Tstar_sph
squarenorm_Tstar_sph !< square of the euclidean norm of Tstar_sph
integer(pInt) :: &
k, l, m, n
associate(prm => param(instance), stt => state(instance))
2019-01-07 11:54:02 +05:30
2018-12-30 17:05:26 +05:30
Tstar_sph = math_spherical33(Tstar)
squarenorm_Tstar_sph = math_mul33xx33(Tstar_sph,Tstar_sph)
2019-01-07 11:54:02 +05:30
norm_Tstar_sph = sqrt(squarenorm_Tstar_sph)
2018-12-30 19:07:31 +05:30
if (prm%dilatation .and. norm_Tstar_sph > 0.0_pReal) then ! no stress or J2 plastitiy --> Li and its derivative are zero
gamma_dot = prm%gdot0 * (sqrt(1.5_pReal) * norm_Tstar_sph /(prm%fTaylor*stt%flowstress(of))) **prm%n
2016-04-22 15:35:06 +05:30
2018-12-30 17:05:26 +05:30
Li = Tstar_sph/norm_Tstar_sph * gamma_dot/prm%fTaylor
2016-04-22 15:35:06 +05:30
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) &
2018-12-30 17:05:26 +05:30
dLi_dTstar(k,l,m,n) = (prm%n-1.0_pReal) * Tstar_sph(k,l)*Tstar_sph(m,n) / squarenorm_Tstar_sph
2016-04-22 15:35:06 +05:30
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt) &
2018-12-30 17:05:26 +05:30
dLi_dTstar(k,l,k,l) = dLi_dTstar(k,l,k,l) + 1.0_pReal
2016-04-22 15:35:06 +05:30
2018-12-30 17:05:26 +05:30
dLi_dTstar = gamma_dot / prm%fTaylor * dLi_dTstar / norm_Tstar_sph
else
2018-12-30 19:07:31 +05:30
Li = 0.0_pReal
dLi_dTstar = 0.0_pReal
endif
2019-01-06 04:11:13 +05:30
2018-12-30 16:03:27 +05:30
end associate
2019-01-06 04:11:13 +05:30
end subroutine plastic_isotropic_LiAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief calculates the rate of change of microstructure
!--------------------------------------------------------------------------------------------------
2018-12-30 19:44:43 +05:30
subroutine plastic_isotropic_dotState(Mp,instance,of)
2016-05-29 14:15:03 +05:30
use prec, only: &
dEq0
use math, only: &
2018-12-30 19:07:31 +05:30
math_mul33xx33, &
math_deviatoric33
2018-12-30 19:44:43 +05:30
implicit none
2018-12-30 19:44:43 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer(pInt), intent(in) :: &
instance, &
of
real(pReal) :: &
gamma_dot, & !< strainrate
hardening, & !< hardening coefficient
saturation, & !< saturation flowstress
2018-12-30 20:39:51 +05:30
norm_Mp !< norm of the (deviatoric) Mandel stress
2018-12-30 19:44:43 +05:30
2018-12-30 19:07:31 +05:30
associate(prm => param(instance), stt => state(instance), dot => dotState(instance))
2019-01-07 11:54:02 +05:30
if (prm%dilatation) then
2018-12-30 19:07:31 +05:30
norm_Mp = sqrt(math_mul33xx33(Mp,Mp))
else
2018-12-30 19:07:31 +05:30
norm_Mp = sqrt(math_mul33xx33(math_deviatoric33(Mp),math_deviatoric33(Mp)))
endif
2019-01-07 11:54:02 +05:30
2018-12-30 19:07:31 +05:30
gamma_dot = prm%gdot0 * (sqrt(1.5_pReal) * norm_Mp /(prm%fTaylor*stt%flowstress(of))) **prm%n
2019-01-07 11:54:02 +05:30
if (abs(gamma_dot) > 1e-12_pReal) then
if (dEq0(prm%tausat_SinhFitA)) then
saturation = prm%tausat
else
saturation = prm%tausat &
2018-12-30 19:07:31 +05:30
+ asinh( (gamma_dot / prm%tausat_SinhFitA)**(1.0_pReal / prm%tausat_SinhFitD) &
)**(1.0_pReal / prm%tausat_SinhFitC) &
2018-12-30 19:07:31 +05:30
/ prm%tausat_SinhFitB * (gamma_dot / prm%gdot0)**(1.0_pReal / prm%n)
endif
hardening = ( prm%h0 + prm%h0_slopeLnRate * log(gamma_dot) ) &
2018-12-30 19:07:31 +05:30
* abs( 1.0_pReal - stt%flowstress(of)/saturation )**prm%a &
* sign(1.0_pReal, 1.0_pReal - stt%flowstress(of)/saturation)
else
hardening = 0.0_pReal
endif
2018-12-30 19:07:31 +05:30
dot%flowstress (of) = hardening * gamma_dot
dot%accumulatedShear(of) = gamma_dot
2019-01-07 11:54:02 +05:30
2018-12-30 16:03:27 +05:30
end associate
end subroutine plastic_isotropic_dotState
2018-12-30 19:07:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
2018-12-30 19:44:43 +05:30
function plastic_isotropic_postResults(Mp,instance,of) result(postResults)
use math, only: &
math_mul33xx33, &
math_deviatoric33
implicit none
2018-12-30 19:44:43 +05:30
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer(pInt), intent(in) :: &
instance, &
of
2018-12-30 19:44:43 +05:30
real(pReal), dimension(sum(plastic_isotropic_sizePostResult(:,instance))) :: &
postResults
real(pReal) :: &
norm_Mp !< norm of the Mandel stress
integer(pInt) :: &
2018-12-30 19:44:43 +05:30
o,c
2018-12-30 19:44:43 +05:30
associate(prm => param(instance), stt => state(instance))
2019-01-07 11:54: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-01-07 11:54:02 +05:30
c = 0_pInt
2018-12-30 16:03:27 +05:30
outputsLoop: do o = 1_pInt,size(prm%outputID)
select case(prm%outputID(o))
2019-01-06 04:11:13 +05:30
case (flowstress_ID)
2018-12-30 19:44:43 +05:30
postResults(c+1_pInt) = stt%flowstress(of)
c = c + 1_pInt
case (strainrate_ID)
2018-12-30 19:44:43 +05:30
postResults(c+1_pInt) = prm%gdot0 &
* (sqrt(1.5_pReal) * norm_Mp /(prm%fTaylor * stt%flowstress(of)))**prm%n
c = c + 1_pInt
2019-01-06 04:11:13 +05:30
end select
enddo outputsLoop
2019-01-07 11:54:02 +05:30
2018-12-30 16:03:27 +05:30
end associate
end function plastic_isotropic_postResults
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine plastic_isotropic_results(instance,group)
#if defined(PETSc) || defined(DAMASKHDF5)
use results
implicit none
integer, intent(in) :: instance
character(len=*) :: group
integer :: o
associate(prm => param(instance), stt => state(instance))
outputsLoop: do o = 1_pInt,size(prm%outputID)
select case(prm%outputID(o))
end select
enddo outputsLoop
end associate
#else
integer, intent(in) :: instance
character(len=*) :: group
#endif
end subroutine plastic_isotropic_results
end module plastic_isotropic