DAMASK_EICMD/src/phase_damage_anisobrittle.f90

213 lines
7.9 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Luv Sharma, Max-Planck-Institut für Eisenforschung GmbH
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine incorporating anisotropic brittle damage source mechanism
!> @details to be done
!--------------------------------------------------------------------------------------------------
submodule (phase:damage) anisobrittle
2020-02-29 19:04:19 +05:30
type :: tParameters !< container type for internal constitutive parameters
real(pREAL) :: &
dot_o_0, & !< opening rate of cleavage planes
p !< damage rate sensitivity
real(pREAL), dimension(:), allocatable :: &
2021-01-19 15:09:16 +05:30
s_crit, & !< critical displacement
2020-09-23 04:38:13 +05:30
g_crit !< critical load
real(pREAL), dimension(:,:,:,:), allocatable :: &
cleavage_systems
integer :: &
2020-08-15 19:32:10 +05:30
sum_N_cl !< total number of cleavage planes
2023-06-04 10:47:38 +05:30
character(len=pSTRLEN), allocatable, dimension(:) :: &
output
end type tParameters
type(tParameters), dimension(:), allocatable :: param !< containers of constitutive parameters (len Ninstances)
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2021-02-13 17:14:47 +05:30
module function anisobrittle_init() result(mySources)
2021-02-13 17:14:47 +05:30
logical, dimension(:), allocatable :: mySources
2020-08-15 19:32:10 +05:30
type(tDict), pointer :: &
2020-08-15 19:32:10 +05:30
phases, &
phase, &
src
2021-07-17 15:20:21 +05:30
integer :: Nmembers,ph
2020-03-17 04:01:43 +05:30
integer, dimension(:), allocatable :: N_cl
character(len=:), allocatable :: &
refs, &
extmsg
2021-02-13 17:14:47 +05:30
mySources = source_active('anisobrittle')
if (count(mySources) == 0) return
2021-02-13 17:37:35 +05:30
print'(/,1x,a)', '<<<+- phase:damage:anisobrittle init -+>>>'
2023-08-29 07:11:48 +05:30
print'(/,1x,a,1x,i0)', '# phases:',count(mySources); flush(IO_STDOUT)
2021-02-13 17:37:35 +05:30
phases => config_material%get_dict('phase')
2021-02-13 14:41:39 +05:30
allocate(param(phases%length))
extmsg = ''
2020-08-15 19:32:10 +05:30
2021-07-17 15:20:21 +05:30
do ph = 1, phases%length
if (mySources(ph)) then
phase => phases%get_dict(ph)
src => phase%get_dict('damage')
2021-02-13 17:14:47 +05:30
2021-07-17 15:20:21 +05:30
associate(prm => param(ph))
2021-01-19 15:09:16 +05:30
2023-08-29 07:11:48 +05:30
print'(/,1x,a,1x,i0,a)', 'phase',ph,': '//phases%key(ph)
refs = config_listReferences(src,indent=3)
if (len(refs) > 0) print'(/,1x,a)', refs
2021-03-11 22:30:07 +05:30
N_cl = src%get_as1dInt('N_cl',defaultVal=emptyIntArray)
2020-08-15 19:32:10 +05:30
prm%sum_N_cl = sum(abs(N_cl))
2021-01-19 15:09:16 +05:30
prm%p = src%get_asReal('p')
prm%dot_o_0 = src%get_asReal('dot_o_0')
2021-01-19 15:09:16 +05:30
prm%s_crit = src%get_as1dReal('s_crit',requiredSize=size(N_cl))
prm%g_crit = src%get_as1dReal('g_crit',requiredSize=size(N_cl))
2021-01-19 15:09:16 +05:30
2023-07-15 23:58:57 +05:30
prm%cleavage_systems = crystal_SchmidMatrix_cleavage(N_cl,phase_lattice(ph),phase_cOverA(ph))
2021-01-19 15:09:16 +05:30
2020-08-15 19:32:10 +05:30
! expand: family => system
2020-09-23 04:38:13 +05:30
prm%s_crit = math_expand(prm%s_crit,N_cl)
prm%g_crit = math_expand(prm%g_crit,N_cl)
2020-08-15 19:32:10 +05:30
#if defined (__GFORTRAN__)
2023-06-04 10:47:38 +05:30
prm%output = output_as1dStr(src)
2020-08-15 19:32:10 +05:30
#else
2023-06-04 10:47:38 +05:30
prm%output = src%get_as1dStr('output',defaultVal=emptyStrArray)
2020-08-15 19:32:10 +05:30
#endif
2021-01-19 15:09:16 +05:30
! sanity checks
if (prm%p <= 0.0_pREAL) extmsg = trim(extmsg)//' p'
if (prm%dot_o_0 <= 0.0_pREAL) extmsg = trim(extmsg)//' dot_o_0'
if (any(prm%g_crit < 0.0_pREAL)) extmsg = trim(extmsg)//' g_crit'
if (any(prm%s_crit < 0.0_pREAL)) extmsg = trim(extmsg)//' s_crit'
2020-08-15 19:32:10 +05:30
2023-01-23 13:01:59 +05:30
Nmembers = count(material_ID_phase==ph)
2021-07-17 15:20:21 +05:30
call phase_allocateState(damageState(ph),Nmembers,1,1,0)
damageState(ph)%atol = src%get_asReal('atol_phi',defaultVal=1.0e-9_pREAL)
if (any(damageState(ph)%atol < 0.0_pREAL)) extmsg = trim(extmsg)//' atol_phi'
2020-08-15 19:32:10 +05:30
2021-07-17 02:11:38 +05:30
end associate
2020-03-15 00:33:57 +05:30
2021-07-17 02:11:38 +05:30
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(damage_anisoBrittle)')
end if
2021-02-13 17:14:47 +05:30
end do
2021-01-26 04:08:32 +05:30
end function anisobrittle_init
!--------------------------------------------------------------------------------------------------
2021-06-25 16:56:08 +05:30
!> @brief
!--------------------------------------------------------------------------------------------------
module subroutine anisobrittle_dotState(M_i, ph,en)
integer, intent(in) :: &
ph,en
real(pREAL), intent(in), dimension(3,3) :: &
M_i
2020-02-29 11:55:36 +05:30
integer :: &
2023-03-10 04:39:21 +05:30
a, i
real(pREAL) :: &
2023-03-10 04:39:21 +05:30
traction, traction_crit
2021-02-13 14:41:39 +05:30
associate(prm => param(ph))
damageState(ph)%dotState(1,en) = 0.0_pREAL
2023-03-10 04:23:47 +05:30
do a = 1, prm%sum_N_cl
traction_crit = damage_phi(ph,en)**2 * prm%g_crit(a)
2023-03-10 04:39:21 +05:30
do i = 1,3
traction = math_tensordot(M_i,prm%cleavage_systems(1:3,1:3,i,a))
2023-03-10 04:39:21 +05:30
damageState(ph)%dotState(1,en) = damageState(ph)%dotState(1,en) &
+ prm%dot_o_0 / prm%s_crit(a) &
* (max(0.0_pREAL, abs(traction) - traction_crit)/traction_crit)**prm%p
2023-03-10 04:39:21 +05:30
end do
end do
2020-02-29 11:55:36 +05:30
end associate
2021-01-26 04:08:32 +05:30
end subroutine anisobrittle_dotState
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
!> @brief Write results to HDF5 output file.
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
module subroutine anisobrittle_result(phase,group)
integer, intent(in) :: phase
2020-02-26 23:07:17 +05:30
character(len=*), intent(in) :: group
integer :: o
2021-02-13 15:31:08 +05:30
2021-02-13 14:41:39 +05:30
associate(prm => param(phase), stt => damageState(phase)%state)
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
case ('Psi_D')
2023-03-22 03:11:52 +05:30
call result_writeDataset(stt,group,trim(prm%output(o)),'damage energy density','J/m³')
2021-02-13 14:41:39 +05:30
end select
end do outputsLoop
end associate
2023-01-19 22:07:45 +05:30
end subroutine anisobrittle_result
!--------------------------------------------------------------------------------------------------
!> @brief contains the constitutive equation for calculating the velocity gradient
!--------------------------------------------------------------------------------------------------
module subroutine damage_anisobrittle_LiAndItsTangent(L_i, dL_i_dM_i, M_i, ph,en)
integer, intent(in) :: &
ph,en
real(pREAL), intent(in), dimension(3,3) :: &
M_i
real(pREAL), intent(out), dimension(3,3) :: &
L_i !< damage velocity gradient
real(pREAL), intent(out), dimension(3,3,3,3) :: &
dL_i_dM_i !< derivative of L_i with respect to M_i
integer :: &
2023-03-10 04:39:21 +05:30
a, k, l, m, n, i
real(pREAL) :: &
2023-03-10 04:23:47 +05:30
traction, traction_crit, &
udot, dudot_dt
L_i = 0.0_pREAL
dL_i_dM_i = 0.0_pREAL
associate(prm => param(ph))
2023-03-10 04:23:47 +05:30
do a = 1,prm%sum_N_cl
traction_crit = damage_phi(ph,en)**2 * prm%g_crit(a)
2023-03-10 04:23:47 +05:30
2023-03-10 04:39:21 +05:30
do i = 1, 3
traction = math_tensordot(M_i,prm%cleavage_systems(1:3,1:3,i,a))
2023-03-10 04:39:21 +05:30
if (abs(traction) > traction_crit + tol_math_check) then
udot = sign(1.0_pREAL,traction)* prm%dot_o_0 * ((abs(traction) - traction_crit)/traction_crit)**prm%p
L_i = L_i + udot*prm%cleavage_systems(1:3,1:3,i,a)
dudot_dt = sign(1.0_pREAL,traction)*udot*prm%p / (abs(traction) - traction_crit)
2023-03-10 04:39:21 +05:30
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dL_i_dM_i(k,l,m,n) = dL_i_dM_i(k,l,m,n) &
+ dudot_dt*prm%cleavage_systems(k,l,i,a) * prm%cleavage_systems(m,n,i,a)
2023-03-10 04:39:21 +05:30
end if
end do
2023-03-10 04:23:47 +05:30
end do
end associate
2021-04-06 15:48:48 +05:30
end subroutine damage_anisobrittle_LiAndItsTangent
end submodule anisobrittle