DAMASK_EICMD/src/damage_nonlocal.f90

240 lines
8.9 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine for non-locally evolving damage field
!--------------------------------------------------------------------------------------------------
module damage_nonlocal
use prec
use material
use config
2019-12-21 15:25:11 +05:30
use numerics
use crystallite
use lattice
use source_damage_isoBrittle
use source_damage_isoDuctile
use source_damage_anisoBrittle
use source_damage_anisoDuctile
2019-12-10 11:44:39 +05:30
use results
implicit none
private
enum, bind(c)
2019-06-11 13:18:07 +05:30
enumerator :: &
undefined_ID, &
damage_ID
end enum
type :: tParameters
2019-06-11 13:18:07 +05:30
integer(kind(undefined_ID)), dimension(:), allocatable :: &
outputID
end type tParameters
2019-06-11 13:18:07 +05:30
type(tparameters), dimension(:), allocatable :: &
param
public :: &
damage_nonlocal_init, &
damage_nonlocal_getSourceAndItsTangent, &
damage_nonlocal_getDiffusion33, &
damage_nonlocal_getMobility, &
damage_nonlocal_putNonLocalDamage, &
2019-12-10 11:44:39 +05:30
damage_nonlocal_Results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine damage_nonlocal_init
2019-12-21 15:25:11 +05:30
integer :: maxNinstance,o,NofMyHomog,h
2019-12-21 17:07:02 +05:30
character(len=pStringLen), dimension(:), allocatable :: outputs
2019-12-21 14:50:50 +05:30
write(6,'(/,a)') ' <<<+- damage_'//DAMAGE_nonlocal_label//' init -+>>>'; flush(6)
maxNinstance = count(damage_type == DAMAGE_nonlocal_ID)
if (maxNinstance == 0) return
allocate(param(maxNinstance))
do h = 1, size(damage_type)
if (damage_type(h) /= DAMAGE_NONLOCAL_ID) cycle
2019-12-21 15:25:11 +05:30
associate(prm => param(damage_typeInstance(h)),config => config_homogenization(h))
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
allocate(prm%outputID(0))
2019-12-21 15:25:11 +05:30
do o=1, size(outputs)
select case(outputs(o))
case ('damage')
prm%outputID = [prm%outputID, damage_ID]
end select
enddo
2019-12-21 15:25:11 +05:30
NofMyHomog = count(material_homogenizationAt == h)
damageState(h)%sizeState = 1
allocate(damageState(h)%state0 (1,NofMyHomog), source=damage_initialPhi(h))
allocate(damageState(h)%subState0(1,NofMyHomog), source=damage_initialPhi(h))
allocate(damageState(h)%state (1,NofMyHomog), source=damage_initialPhi(h))
2019-12-21 15:25:11 +05:30
nullify(damageMapping(h)%p)
damageMapping(h)%p => material_homogenizationMemberAt
2019-12-21 15:25:11 +05:30
deallocate(damage(h)%p)
damage(h)%p => damageState(h)%state(1,:)
end associate
enddo
2019-12-21 15:25:11 +05:30
end subroutine damage_nonlocal_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates homogenized damage driving forces
!--------------------------------------------------------------------------------------------------
subroutine damage_nonlocal_getSourceAndItsTangent(phiDot, dPhiDot_dPhi, phi, ip, el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: &
phi
integer :: &
phase, &
grain, &
source, &
constituent
real(pReal) :: &
phiDot, dPhiDot_dPhi, localphiDot, dLocalphiDot_dPhi
phiDot = 0.0_pReal
dPhiDot_dPhi = 0.0_pReal
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
phase = material_phaseAt(grain,el)
constituent = material_phasememberAt(grain,ip,el)
do source = 1, phase_Nsources(phase)
select case(phase_source(source,phase))
case (SOURCE_damage_isoBrittle_ID)
call source_damage_isobrittle_getRateAndItsTangent (localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
case (SOURCE_damage_isoDuctile_ID)
call source_damage_isoductile_getRateAndItsTangent (localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
case (SOURCE_damage_anisoBrittle_ID)
call source_damage_anisobrittle_getRateAndItsTangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
case (SOURCE_damage_anisoDuctile_ID)
call source_damage_anisoductile_getRateAndItsTangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
case default
localphiDot = 0.0_pReal
dLocalphiDot_dPhi = 0.0_pReal
end select
phiDot = phiDot + localphiDot
dPhiDot_dPhi = dPhiDot_dPhi + dLocalphiDot_dPhi
enddo
enddo
phiDot = phiDot/real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
dPhiDot_dPhi = dPhiDot_dPhi/real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
end subroutine damage_nonlocal_getSourceAndItsTangent
2019-05-17 02:44:47 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized non local damage diffusion tensor in reference configuration
!--------------------------------------------------------------------------------------------------
function damage_nonlocal_getDiffusion33(ip,el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
damage_nonlocal_getDiffusion33
integer :: &
homog, &
grain
homog = material_homogenizationAt(el)
damage_nonlocal_getDiffusion33 = 0.0_pReal
do grain = 1, homogenization_Ngrains(homog)
damage_nonlocal_getDiffusion33 = damage_nonlocal_getDiffusion33 + &
2019-06-15 17:27:24 +05:30
crystallite_push33ToRef(grain,ip,el,lattice_DamageDiffusion33(1:3,1:3,material_phaseAt(grain,el)))
enddo
damage_nonlocal_getDiffusion33 = &
charLength**2*damage_nonlocal_getDiffusion33/real(homogenization_Ngrains(homog),pReal)
end function damage_nonlocal_getDiffusion33
2019-05-17 02:44:47 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Returns homogenized nonlocal damage mobility
!--------------------------------------------------------------------------------------------------
real(pReal) function damage_nonlocal_getMobility(ip,el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
integer :: &
ipc
damage_nonlocal_getMobility = 0.0_pReal
2019-06-07 00:44:37 +05:30
do ipc = 1, homogenization_Ngrains(material_homogenizationAt(el))
2019-06-15 17:27:24 +05:30
damage_nonlocal_getMobility = damage_nonlocal_getMobility + lattice_DamageMobility(material_phaseAt(ipc,el))
enddo
damage_nonlocal_getMobility = damage_nonlocal_getMobility/&
2019-06-07 00:44:37 +05:30
real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
end function damage_nonlocal_getMobility
2019-05-17 02:44:47 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief updated nonlocal damage field with solution from damage phase field PDE
!--------------------------------------------------------------------------------------------------
subroutine damage_nonlocal_putNonLocalDamage(phi,ip,el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: &
phi
integer :: &
homog, &
offset
homog = material_homogenizationAt(el)
offset = damageMapping(homog)%p(ip,el)
damage(homog)%p(offset) = phi
end subroutine damage_nonlocal_putNonLocalDamage
2019-05-17 02:44:47 +05:30
2019-12-10 11:44:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
2019-12-10 21:55:51 +05:30
subroutine damage_nonlocal_results(homog,group)
2019-12-10 11:44:39 +05:30
integer, intent(in) :: homog
character(len=*), intent(in) :: group
2019-12-21 15:06:42 +05:30
integer :: o
2019-12-10 11:44:39 +05:30
2019-12-21 15:06:42 +05:30
associate(prm => param(damage_typeInstance(homog)))
2019-12-10 11:44:39 +05:30
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
case (damage_ID)
call results_writeDataset(group,damage(homog)%p,'phi',&
'damage indicator','-')
end select
enddo outputsLoop
end associate
end subroutine damage_nonlocal_results
end module damage_nonlocal