DAMASK_EICMD/src/phase_damage_anisobrittle.f90

197 lines
7.3 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:damagee) 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
class(tNode), pointer :: &
phases, &
phase, &
sources, &
src
2021-02-13 17:37:35 +05:30
integer :: Nconstituents,p
2020-03-17 04:01:43 +05:30
integer, dimension(:), allocatable :: N_cl
2020-02-29 14:50:38 +05:30
character(len=pStringLen) :: extmsg = ''
2021-02-13 17:14:47 +05:30
mySources = source_active('anisobrittle')
2021-02-13 17:37:35 +05:30
if(count(mySources) == 0) return
print'(/,a)', ' <<<+- phase:damage:anisobrittle init -+>>>'
print'(a,i0)', ' # phases: ',count(mySources); flush(IO_STDOUT)
phases => config_material%get('phase')
2021-02-13 14:41:39 +05:30
allocate(param(phases%length))
2020-08-15 19:32:10 +05:30
do p = 1, phases%length
2021-02-13 17:14:47 +05:30
if(mySources(p)) then
2021-01-19 15:09:16 +05:30
phase => phases%get(p)
sources => phase%get('damage')
2021-02-13 17:14:47 +05:30
2021-02-13 14:41:39 +05:30
associate(prm => param(p))
2021-02-13 17:14:47 +05:30
src => sources%get(1)
2021-01-19 15:09:16 +05:30
2020-08-15 19:32:10 +05:30
N_cl = src%get_asInts('N_cl',defaultVal=emptyIntArray)
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
2020-09-23 04:38:13 +05:30
prm%s_crit = src%get_asFloats('s_crit', requiredSize=size(N_cl))
prm%g_crit = src%get_asFloats('g_crit', requiredSize=size(N_cl))
2021-01-19 15:09:16 +05:30
2020-08-15 19:32:10 +05:30
prm%cleavage_systems = lattice_SchmidMatrix_cleavage(N_cl,phase%get_asString('lattice'),&
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
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__)
prm%output = output_asStrings(src)
#else
prm%output = src%get_asStrings('output',defaultVal=emptyStringArray)
#endif
2021-01-19 15:09:16 +05:30
2020-08-15 19:32:10 +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
Nconstituents = count(material_phaseAt==p) * discretization_nIPs
call phase_allocateState(damageState(p),Nconstituents,1,1,0)
damageState(p)%atol = src%get_asFloat('anisobrittle_atol',defaultVal=1.0e-3_pReal)
if(any(damageState(p)%atol < 0.0_pReal)) extmsg = trim(extmsg)//' anisobrittle_atol'
2020-08-15 19:32:10 +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)//'(damage_anisoBrittle)')
endif
2021-02-13 17:14:47 +05:30
2020-08-15 19:32:10 +05:30
enddo
2021-01-26 04:08:32 +05:30
end function anisobrittle_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
2021-02-13 15:31:08 +05:30
module subroutine anisobrittle_dotState(S, ph,me)
integer, intent(in) :: &
2021-02-13 15:31:08 +05:30
ph,me
real(pReal), intent(in), dimension(3,3) :: &
S
2020-02-29 11:55:36 +05:30
integer :: &
sourceOffset, &
damageOffset, &
homog, &
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,me) = 0.0_pReal
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
2021-02-13 15:31:08 +05:30
traction_crit = prm%g_crit(i)*damage_phi(ph,me)**2.0_pReal
2021-02-13 14:41:39 +05:30
damageState(ph)%dotState(1,me) = damageState(ph)%dotState(1,me) &
+ 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)
enddo
2020-02-29 11:55:36 +05:30
end associate
2021-01-26 04:08:32 +05:30
end subroutine anisobrittle_dotState
!--------------------------------------------------------------------------------------------------
!> @brief returns local part of nonlocal damage driving force
!--------------------------------------------------------------------------------------------------
2021-02-13 15:31:08 +05:30
module subroutine anisobrittle_getRateAndItsTangent(localphiDot, dLocalphiDot_dPhi, phi, ph, me)
integer, intent(in) :: &
2021-02-13 15:31:08 +05:30
ph, &
me
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
localphiDot, &
dLocalphiDot_dPhi
2020-02-29 11:55:36 +05:30
2020-02-26 23:07:17 +05:30
2021-02-13 15:31:08 +05:30
dLocalphiDot_dPhi = -damageState(ph)%state(1,me)
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
localphiDot = 1.0_pReal &
+ dLocalphiDot_dPhi*phi
2021-02-13 15:31:08 +05:30
end subroutine anisobrittle_getRateAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
2021-01-26 04:08:32 +05:30
module subroutine anisobrittle_results(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')
call results_writeDataset(group,stt,trim(prm%output(o)),'driving force','J/m³')
end select
enddo outputsLoop
end associate
2021-01-26 04:08:32 +05:30
end subroutine anisobrittle_results
end submodule anisobrittle