DAMASK_EICMD/src/source_damage_anisoBrittle.f90

228 lines
9.0 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
!--------------------------------------------------------------------------------------------------
module source_damage_anisoBrittle
use prec
use debug
use IO
use math
use discretization
2020-02-29 11:55:36 +05:30
use material
use config
use lattice
use results
implicit none
private
integer, dimension(:), allocatable :: &
source_damage_anisoBrittle_offset, & !< which source is my current source mechanism?
source_damage_anisoBrittle_instance !< instance of source mechanism
2020-02-26 23:07:17 +05:30
2020-02-29 19:04:19 +05:30
type :: tParameters !< container type for internal constitutive parameters
real(pReal) :: &
sdot_0, &
N
real(pReal), dimension(:), allocatable :: &
critDisp, &
critLoad
real(pReal), dimension(:,:,:,:), allocatable :: &
cleavage_systems
integer :: &
totalNcleavage
integer, dimension(:), allocatable :: &
Ncleavage
character(len=pStringLen), allocatable, dimension(:) :: &
output
end type tParameters
type(tParameters), dimension(:), allocatable :: param !< containers of constitutive parameters (len Ninstance)
public :: &
source_damage_anisoBrittle_init, &
source_damage_anisoBrittle_dotState, &
source_damage_anisobrittle_getRateAndItsTangent, &
source_damage_anisoBrittle_results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2019-02-13 14:27:12 +05:30
subroutine source_damage_anisoBrittle_init
2020-02-29 12:33:06 +05:30
integer :: Ninstance,sourceOffset,NofMyPhase,p
2020-02-29 14:50:38 +05:30
character(len=pStringLen) :: extmsg = ''
2020-02-26 23:07:17 +05:30
write(6,'(/,a)') ' <<<+- source_'//SOURCE_DAMAGE_ANISOBRITTLE_LABEL//' init -+>>>'; flush(6)
2020-02-26 23:07:17 +05:30
Ninstance = count(phase_source == SOURCE_DAMAGE_ANISOBRITTLE_ID)
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2020-02-26 23:07:17 +05:30
allocate(source_damage_anisoBrittle_offset (size(config_phase)), source=0)
allocate(source_damage_anisoBrittle_instance(size(config_phase)), source=0)
allocate(param(Ninstance))
2020-02-26 23:07:17 +05:30
do p = 1, size(config_phase)
source_damage_anisoBrittle_instance(p) = count(phase_source(:,1:p) == SOURCE_DAMAGE_ANISOBRITTLE_ID)
2020-02-29 12:33:06 +05:30
do sourceOffset = 1, phase_Nsources(p)
if (phase_source(sourceOffset,p) == SOURCE_DAMAGE_ANISOBRITTLE_ID) then
2020-02-29 12:33:06 +05:30
source_damage_anisoBrittle_offset(p) = sourceOffset
exit
endif
2020-02-26 23:07:17 +05:30
enddo
if (all(phase_source(:,p) /= SOURCE_DAMAGE_ANISOBRITTLE_ID)) cycle
associate(prm => param(source_damage_anisoBrittle_instance(p)), &
config => config_phase(p))
2020-02-26 23:07:17 +05:30
2020-03-15 00:33:57 +05:30
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
2020-02-29 11:55:36 +05:30
prm%Ncleavage = config%getInts('ncleavage',defaultVal=emptyIntArray)
prm%totalNcleavage = sum(prm%Ncleavage)
prm%N = config%getFloat('anisobrittle_ratesensitivity')
prm%sdot_0 = config%getFloat('anisobrittle_sdot0')
2020-02-29 11:55:36 +05:30
2020-02-28 14:55:07 +05:30
prm%critDisp = config%getFloats('anisobrittle_criticaldisplacement',requiredSize=size(prm%Ncleavage))
prm%critLoad = config%getFloats('anisobrittle_criticalload', requiredSize=size(prm%Ncleavage))
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
prm%cleavage_systems = lattice_SchmidMatrix_cleavage(prm%Ncleavage,config%getString('lattice_structure'),&
config%getFloat('c/a',defaultVal=0.0_pReal))
2020-02-26 23:07:17 +05:30
! expand: family => system
2020-02-29 11:55:36 +05:30
prm%critDisp = math_expand(prm%critDisp, prm%Ncleavage)
prm%critLoad = math_expand(prm%critLoad, prm%Ncleavage)
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
! sanity checks
2020-02-29 14:50:38 +05:30
if (prm%N <= 0.0_pReal) extmsg = trim(extmsg)//' anisobrittle_n'
2020-02-29 11:55:36 +05:30
if (prm%sdot_0 <= 0.0_pReal) extmsg = trim(extmsg)//' anisobrittle_sdot0'
2020-02-29 14:50:38 +05:30
if (any(prm%critLoad < 0.0_pReal)) extmsg = trim(extmsg)//' anisobrittle_critLoad'
if (any(prm%critDisp < 0.0_pReal)) extmsg = trim(extmsg)//' anisobrittle_critDisp'
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
NofMyPhase = count(material_phaseAt==p) * discretization_nIP
2020-02-26 23:07:17 +05:30
call material_allocateSourceState(p,sourceOffset,NofMyPhase,1,1,0)
2020-03-15 00:33:57 +05:30
sourceState(p)%p(sourceOffset)%atolState = config%getFloat('anisobrittle_atol',defaultVal=1.0e-3_pReal)
if(any(sourceState(p)%p(sourceOffset)%atolState <= 0.0_pReal)) &
extmsg = trim(extmsg)//' anisobrittle_atol'
2020-02-28 14:55:07 +05:30
end associate
2020-03-15 00:33:57 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'('//SOURCE_DAMAGE_ANISOBRITTLE_LABEL//')')
enddo
end subroutine source_damage_anisoBrittle_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisoBrittle_dotState(S, ipc, ip, el)
integer, intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
real(pReal), intent(in), dimension(3,3) :: &
S
2020-02-29 11:55:36 +05:30
integer :: &
phase, &
constituent, &
sourceOffset, &
damageOffset, &
homog, &
2020-02-29 11:55:36 +05:30
i
real(pReal) :: &
traction_d, traction_t, traction_n, traction_crit
phase = material_phaseAt(ipc,el)
constituent = material_phasememberAt(ipc,ip,el)
sourceOffset = source_damage_anisoBrittle_offset(phase)
homog = material_homogenizationAt(el)
damageOffset = damageMapping(homog)%p(ip,el)
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
associate(prm => param(source_damage_anisoBrittle_instance(phase)))
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) = 0.0_pReal
2020-02-29 11:55:36 +05:30
do i = 1, prm%totalNcleavage
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
traction_d = math_mul33xx33(S,prm%cleavage_systems(1:3,1:3,1,i))
traction_t = math_mul33xx33(S,prm%cleavage_systems(1:3,1:3,2,i))
traction_n = math_mul33xx33(S,prm%cleavage_systems(1:3,1:3,3,i))
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
traction_crit = prm%critLoad(i)*damage(homog)%p(damageOffset)**2.0_pReal
2020-02-29 11:55:36 +05:30
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) &
= sourceState(phase)%p(sourceOffset)%dotState(1,constituent) &
2020-02-29 19:04:19 +05:30
+ prm%sdot_0 / prm%critDisp(i) &
* ((max(0.0_pReal, abs(traction_d) - traction_crit)/traction_crit)**prm%N + &
(max(0.0_pReal, abs(traction_t) - traction_crit)/traction_crit)**prm%N + &
(max(0.0_pReal, abs(traction_n) - traction_crit)/traction_crit)**prm%N)
enddo
2020-02-29 11:55:36 +05:30
end associate
end subroutine source_damage_anisoBrittle_dotState
!--------------------------------------------------------------------------------------------------
!> @brief returns local part of nonlocal damage driving force
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisobrittle_getRateAndItsTangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
integer, intent(in) :: &
phase, &
constituent
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
localphiDot, &
dLocalphiDot_dPhi
2020-02-29 11:55:36 +05:30
integer :: &
sourceOffset
sourceOffset = source_damage_anisoBrittle_offset(phase)
2020-02-26 23:07:17 +05:30
dLocalphiDot_dPhi = -sourceState(phase)%p(sourceOffset)%state(1,constituent)
2020-02-26 23:07:17 +05:30
2020-02-29 11:55:36 +05:30
localphiDot = 1.0_pReal &
+ dLocalphiDot_dPhi*phi
end subroutine source_damage_anisoBrittle_getRateAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisoBrittle_results(phase,group)
integer, intent(in) :: phase
2020-02-26 23:07:17 +05:30
character(len=*), intent(in) :: group
integer :: o
associate(prm => param(source_damage_anisoBrittle_instance(phase)), &
stt => sourceState(phase)%p(source_damage_anisoBrittle_offset(phase))%state)
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
case ('anisobrittle_drivingforce')
call results_writeDataset(group,stt,'tbd','driving force','tbd')
end select
enddo outputsLoop
end associate
end subroutine source_damage_anisoBrittle_results
end module source_damage_anisoBrittle