DAMASK_EICMD/src/phase_damage_anisobrittle.f90

228 lines
9.2 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) :: &
2020-09-23 04:38:13 +05:30
dot_o, & !< opening rate of cleavage planes
q !< 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
character(len=pStringLen), 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 :: 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 -+>>>'
print'(/,a,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
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
2020-09-23 04:38:13 +05:30
prm%q = src%get_asFloat('q')
prm%dot_o = src%get_asFloat('dot_o')
2021-01-19 15:09:16 +05:30
2021-03-11 22:30:07 +05:30
prm%s_crit = src%get_as1dFloat('s_crit', requiredSize=size(N_cl))
prm%g_crit = src%get_as1dFloat('g_crit', requiredSize=size(N_cl))
2021-01-19 15:09:16 +05:30
2021-07-21 19:16:38 +05:30
prm%cleavage_systems = lattice_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__)
2021-03-11 22:30:07 +05:30
prm%output = output_as1dString(src)
2020-08-15 19:32:10 +05:30
#else
2021-03-11 22:30:07 +05:30
prm%output = src%get_as1dString('output',defaultVal=emptyStringArray)
2020-08-15 19:32:10 +05:30
#endif
2021-01-19 15:09:16 +05:30
! sanity checks
2020-09-23 04:38:13 +05:30
if (prm%q <= 0.0_pReal) extmsg = trim(extmsg)//' q'
if (prm%dot_o <= 0.0_pReal) extmsg = trim(extmsg)//' dot_o'
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
2021-07-17 15:20:21 +05:30
Nmembers = count(material_phaseID==ph)
call phase_allocateState(damageState(ph),Nmembers,1,1,0)
damageState(ph)%atol = src%get_asFloat('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(S, ph,en)
integer, intent(in) :: &
ph,en
real(pReal), intent(in), dimension(3,3) :: &
S
2020-02-29 11:55:36 +05:30
integer :: &
2020-02-29 11:55:36 +05:30
i
real(pReal) :: &
traction_d, traction_t, traction_n, traction_crit
2021-02-13 14:41:39 +05:30
associate(prm => param(ph))
damageState(ph)%dotState(1,en) = 0.0_pReal
2021-02-13 14:41:39 +05:30
do i = 1, prm%sum_N_cl
traction_d = math_tensordot(S,prm%cleavage_systems(1:3,1:3,1,i))
traction_t = math_tensordot(S,prm%cleavage_systems(1:3,1:3,2,i))
traction_n = math_tensordot(S,prm%cleavage_systems(1:3,1:3,3,i))
2020-02-26 23:07:17 +05:30
traction_crit = prm%g_crit(i)*damage_phi(ph,en)**2
damageState(ph)%dotState(1,en) = damageState(ph)%dotState(1,en) &
2021-02-13 14:41:39 +05:30
+ prm%dot_o / prm%s_crit(i) &
* ((max(0.0_pReal, abs(traction_d) - traction_crit)/traction_crit)**prm%q + &
(max(0.0_pReal, abs(traction_t) - traction_crit)/traction_crit)**prm%q + &
(max(0.0_pReal, abs(traction_n) - traction_crit)/traction_crit)**prm%q)
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 ('f_phi')
2023-01-19 22:07:45 +05:30
call result_writeDataset(stt,group,trim(prm%output(o)),'driving force','-')
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(Ld, dLd_dTstar, S, ph,en)
integer, intent(in) :: &
ph,en
real(pReal), intent(in), dimension(3,3) :: &
S
real(pReal), intent(out), dimension(3,3) :: &
Ld !< damage velocity gradient
real(pReal), intent(out), dimension(3,3,3,3) :: &
dLd_dTstar !< derivative of Ld with respect to Tstar (4th-order tensor)
integer :: &
i, k, l, m, n
real(pReal) :: &
traction_d, traction_t, traction_n, traction_crit, &
udotd, dudotd_dt, udott, dudott_dt, udotn, dudotn_dt
Ld = 0.0_pReal
dLd_dTstar = 0.0_pReal
associate(prm => param(ph))
do i = 1,prm%sum_N_cl
traction_crit = prm%g_crit(i)*damage_phi(ph,en)**2
traction_d = math_tensordot(S,prm%cleavage_systems(1:3,1:3,1,i))
if (abs(traction_d) > traction_crit + tol_math_check) then
udotd = sign(1.0_pReal,traction_d)* prm%dot_o * ((abs(traction_d) - traction_crit)/traction_crit)**prm%q
Ld = Ld + udotd*prm%cleavage_systems(1:3,1:3,1,i)
dudotd_dt = sign(1.0_pReal,traction_d)*udotd*prm%q / (abs(traction_d) - traction_crit)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLd_dTstar(k,l,m,n) = dLd_dTstar(k,l,m,n) &
+ dudotd_dt*prm%cleavage_systems(k,l,1,i) * prm%cleavage_systems(m,n,1,i)
end if
traction_t = math_tensordot(S,prm%cleavage_systems(1:3,1:3,2,i))
if (abs(traction_t) > traction_crit + tol_math_check) then
udott = sign(1.0_pReal,traction_t)* prm%dot_o * ((abs(traction_t) - traction_crit)/traction_crit)**prm%q
Ld = Ld + udott*prm%cleavage_systems(1:3,1:3,2,i)
dudott_dt = sign(1.0_pReal,traction_t)*udott*prm%q / (abs(traction_t) - traction_crit)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLd_dTstar(k,l,m,n) = dLd_dTstar(k,l,m,n) &
+ dudott_dt*prm%cleavage_systems(k,l,2,i) * prm%cleavage_systems(m,n,2,i)
end if
traction_n = math_tensordot(S,prm%cleavage_systems(1:3,1:3,3,i))
if (abs(traction_n) > traction_crit + tol_math_check) then
udotn = sign(1.0_pReal,traction_n)* prm%dot_o * ((abs(traction_n) - traction_crit)/traction_crit)**prm%q
Ld = Ld + udotn*prm%cleavage_systems(1:3,1:3,3,i)
dudotn_dt = sign(1.0_pReal,traction_n)*udotn*prm%q / (abs(traction_n) - traction_crit)
forall (k=1:3,l=1:3,m=1:3,n=1:3) &
dLd_dTstar(k,l,m,n) = dLd_dTstar(k,l,m,n) &
+ dudotn_dt*prm%cleavage_systems(k,l,3,i) * prm%cleavage_systems(m,n,3,i)
end if
end do
end associate
2021-04-06 15:48:48 +05:30
end subroutine damage_anisobrittle_LiAndItsTangent
end submodule anisobrittle