DAMASK_EICMD/src/source_damage_anisoBrittle.f90

329 lines
13 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
2019-05-15 02:14:38 +05:30
use prec
implicit none
private
2019-05-15 02:14:38 +05:30
integer, dimension(:), allocatable, public, protected :: &
source_damage_anisoBrittle_offset, & !< which source is my current source mechanism?
source_damage_anisoBrittle_instance !< instance of source mechanism
2019-05-15 02:14:38 +05:30
integer, dimension(:,:), allocatable, target, public :: &
source_damage_anisoBrittle_sizePostResult !< size of each post result output
2019-05-15 02:14:38 +05:30
character(len=64), dimension(:,:), allocatable, target, public :: &
source_damage_anisoBrittle_output !< name of each post result output
2019-05-15 02:14:38 +05:30
integer, dimension(:,:), allocatable, private :: &
source_damage_anisoBrittle_Ncleavage !< number of cleavage systems per family
enum, bind(c)
enumerator :: undefined_ID, &
damage_drivingforce_ID
end enum
type, private :: tParameters !< container type for internal constitutive parameters
real(pReal) :: &
aTol, &
sdot_0, &
N
real(pReal), dimension(:), allocatable :: &
critDisp, &
critLoad
real(pReal), dimension(:,:,:,:), allocatable :: &
cleavage_systems
2019-05-15 02:14:38 +05:30
integer :: &
totalNcleavage
2019-05-15 02:14:38 +05:30
integer, dimension(:), allocatable :: &
Ncleavage
2019-02-13 12:36:22 +05:30
integer(kind(undefined_ID)), allocatable, dimension(:) :: &
outputID !< ID of each post result output
end type tParameters
type(tParameters), dimension(:), allocatable, private :: param !< containers of constitutive parameters (len Ninstance)
public :: &
source_damage_anisoBrittle_init, &
source_damage_anisoBrittle_dotState, &
source_damage_anisobrittle_getRateAndItsTangent, &
source_damage_anisoBrittle_postResults
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
use debug, only: &
debug_level,&
debug_constitutive,&
debug_levelBasic
use IO, only: &
2019-02-13 14:33:28 +05:30
IO_error
2019-02-13 13:46:06 +05:30
use math, only: &
math_expand
use material, only: &
2019-02-13 11:52:37 +05:30
material_allocateSourceState, &
phase_source, &
phase_Nsources, &
phase_Noutput, &
SOURCE_damage_anisoBrittle_label, &
SOURCE_damage_anisoBrittle_ID, &
material_phase, &
sourceState
use config, only: &
config_phase, &
2019-03-09 05:37:26 +05:30
material_Nphase
use lattice, only: &
lattice_SchmidMatrix_cleavage, &
2019-02-13 14:27:12 +05:30
lattice_maxNcleavageFamily
2019-05-15 02:14:38 +05:30
integer :: Ninstance,phase,instance,source,sourceOffset
integer :: NofMyPhase,p ,i
integer, dimension(0), parameter :: emptyIntArray = [integer::]
character(len=65536), dimension(0), parameter :: emptyStringArray = [character(len=65536)::]
2019-02-13 12:36:22 +05:30
integer(kind(undefined_ID)) :: &
outputID
character(len=pStringLen) :: &
extmsg = ''
character(len=65536), dimension(:), allocatable :: &
outputs
2019-02-13 12:36:22 +05:30
write(6,'(/,a)') ' <<<+- source_'//SOURCE_DAMAGE_ANISOBRITTLE_LABEL//' init -+>>>'
2019-05-17 02:44:47 +05:30
Ninstance = count(phase_source == SOURCE_damage_anisoBrittle_ID)
2019-05-15 02:14:38 +05:30
if (Ninstance == 0) return
2019-05-15 02:14:38 +05:30
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-05-15 02:14:38 +05:30
allocate(source_damage_anisoBrittle_offset(material_Nphase), source=0)
allocate(source_damage_anisoBrittle_instance(material_Nphase), source=0)
do phase = 1, material_Nphase
source_damage_anisoBrittle_instance(phase) = count(phase_source(:,1:phase) == source_damage_anisoBrittle_ID)
do source = 1, phase_Nsources(phase)
if (phase_source(source,phase) == source_damage_anisoBrittle_ID) &
source_damage_anisoBrittle_offset(phase) = source
enddo
enddo
2019-02-13 13:35:07 +05:30
2019-05-15 02:14:38 +05:30
allocate(source_damage_anisoBrittle_sizePostResult(maxval(phase_Noutput),Ninstance), source=0)
allocate(source_damage_anisoBrittle_output(maxval(phase_Noutput),Ninstance))
source_damage_anisoBrittle_output = ''
2019-05-15 02:14:38 +05:30
allocate(source_damage_anisoBrittle_Ncleavage(lattice_maxNcleavageFamily,Ninstance), source=0)
allocate(param(Ninstance))
do p=1, size(config_phase)
if (all(phase_source(:,p) /= SOURCE_DAMAGE_ANISOBRITTLE_ID)) cycle
associate(prm => param(source_damage_anisoBrittle_instance(p)), &
config => config_phase(p))
prm%aTol = config%getFloat('anisobrittle_atol',defaultVal = 1.0e-3_pReal)
prm%N = config%getFloat('anisobrittle_ratesensitivity')
prm%sdot_0 = config%getFloat('anisobrittle_sdot0')
! sanity checks
if (prm%aTol < 0.0_pReal) extmsg = trim(extmsg)//' anisobrittle_atol'
if (prm%N <= 0.0_pReal) extmsg = trim(extmsg)//' anisobrittle_ratesensitivity'
if (prm%sdot_0 <= 0.0_pReal) extmsg = trim(extmsg)//' anisobrittle_sdot0'
prm%Ncleavage = config%getInts('ncleavage',defaultVal=emptyIntArray)
2019-02-13 13:46:06 +05:30
prm%critDisp = config%getFloats('anisobrittle_criticaldisplacement',requiredSize=size(prm%Ncleavage))
prm%critLoad = config%getFloats('anisobrittle_criticalload', requiredSize=size(prm%Ncleavage))
prm%cleavage_systems = lattice_SchmidMatrix_cleavage (prm%Ncleavage,config%getString('lattice_structure'),&
config%getFloat('c/a',defaultVal=0.0_pReal))
2019-02-13 13:46:06 +05:30
! expand: family => system
prm%critDisp = math_expand(prm%critDisp, prm%Ncleavage)
prm%critLoad = math_expand(prm%critLoad, prm%Ncleavage)
2019-02-13 14:27:12 +05:30
if (any(prm%critLoad < 0.0_pReal)) extmsg = trim(extmsg)//' anisobrittle_criticalload'
if (any(prm%critDisp < 0.0_pReal)) extmsg = trim(extmsg)//' anisobrittle_criticaldisplacement'
2019-02-13 12:36:22 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if (extmsg /= '') &
2019-05-15 02:14:38 +05:30
call IO_error(211,ext_msg=trim(extmsg)//'('//SOURCE_DAMAGE_ANISOBRITTLE_LABEL//')')
2019-02-13 12:36:22 +05:30
!--------------------------------------------------------------------------------------------------
! output pararameters
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
allocate(prm%outputID(0))
2019-05-15 02:14:38 +05:30
do i=1, size(outputs)
2019-02-13 12:36:22 +05:30
outputID = undefined_ID
select case(outputs(i))
2019-02-13 13:35:07 +05:30
2019-02-13 12:36:22 +05:30
case ('anisobrittle_drivingforce')
2019-05-15 02:14:38 +05:30
source_damage_anisoBrittle_sizePostResult(i,source_damage_anisoBrittle_instance(p)) = 1
2019-02-13 13:35:07 +05:30
source_damage_anisoBrittle_output(i,source_damage_anisoBrittle_instance(p)) = outputs(i)
prm%outputID = [prm%outputID, damage_drivingforce_ID]
2019-02-13 12:36:22 +05:30
end select
enddo
end associate
2019-02-13 13:35:07 +05:30
phase = p
NofMyPhase=count(material_phase==phase)
instance = source_damage_anisoBrittle_instance(phase)
sourceOffset = source_damage_anisoBrittle_offset(phase)
2019-05-15 02:14:38 +05:30
call material_allocateSourceState(phase,sourceOffset,NofMyPhase,1,1,0)
2019-02-13 13:35:07 +05:30
sourceState(phase)%p(sourceOffset)%sizePostResults = sum(source_damage_anisoBrittle_sizePostResult(:,instance))
sourceState(phase)%p(sourceOffset)%aTolState=param(instance)%aTol
2019-02-13 14:27:12 +05:30
source_damage_anisoBrittle_Ncleavage(1:size(param(instance)%Ncleavage),instance) = param(instance)%Ncleavage
enddo
end subroutine source_damage_anisoBrittle_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisoBrittle_dotState(S, ipc, ip, el)
use math, only: &
math_mul33xx33
use material, only: &
phaseAt, phasememberAt, &
sourceState, &
2019-03-10 15:32:32 +05:30
material_homogenizationAt, &
damage, &
damageMapping
use lattice, only: &
lattice_Scleavage, &
lattice_maxNcleavageFamily, &
lattice_NcleavageSystem
2019-05-15 02:14:38 +05:30
integer, intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
real(pReal), intent(in), dimension(3,3) :: &
S
2019-05-15 02:14:38 +05:30
integer :: &
phase, &
constituent, &
instance, &
sourceOffset, &
damageOffset, &
homog, &
2019-02-13 14:27:12 +05:30
f, i, index_myFamily, index
real(pReal) :: &
traction_d, traction_t, traction_n, traction_crit
phase = phaseAt(ipc,ip,el)
constituent = phasememberAt(ipc,ip,el)
instance = source_damage_anisoBrittle_instance(phase)
sourceOffset = source_damage_anisoBrittle_offset(phase)
2019-03-10 15:32:32 +05:30
homog = material_homogenizationAt(el)
damageOffset = damageMapping(homog)%p(ip,el)
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) = 0.0_pReal
2019-02-13 14:27:12 +05:30
2019-05-15 02:14:38 +05:30
index = 1
do f = 1,lattice_maxNcleavageFamily
index_myFamily = sum(lattice_NcleavageSystem(1:f-1,phase)) ! at which index starts my family
do i = 1,source_damage_anisoBrittle_Ncleavage(f,instance) ! process each (active) cleavage system in family
traction_d = math_mul33xx33(S,lattice_Scleavage(1:3,1:3,1,index_myFamily+i,phase))
traction_t = math_mul33xx33(S,lattice_Scleavage(1:3,1:3,2,index_myFamily+i,phase))
traction_n = math_mul33xx33(S,lattice_Scleavage(1:3,1:3,3,index_myFamily+i,phase))
2019-02-13 14:27:12 +05:30
traction_crit = param(instance)%critLoad(index)* &
damage(homog)%p(damageOffset)*damage(homog)%p(damageOffset)
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) = &
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) + &
2019-02-13 12:06:36 +05:30
param(instance)%sdot_0* &
((max(0.0_pReal, abs(traction_d) - traction_crit)/traction_crit)**param(instance)%N + &
(max(0.0_pReal, abs(traction_t) - traction_crit)/traction_crit)**param(instance)%N + &
(max(0.0_pReal, abs(traction_n) - traction_crit)/traction_crit)**param(instance)%N)/ &
2019-02-13 14:27:12 +05:30
param(instance)%critDisp(index)
2019-05-15 02:14:38 +05:30
index = index + 1
enddo
enddo
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)
use material, only: &
sourceState
2019-05-15 02:14:38 +05:30
integer, intent(in) :: &
phase, &
constituent
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
localphiDot, &
dLocalphiDot_dPhi
2019-05-15 02:14:38 +05:30
integer :: &
sourceOffset
sourceOffset = source_damage_anisoBrittle_offset(phase)
2019-02-13 14:46:06 +05:30
localphiDot = 1.0_pReal &
- sourceState(phase)%p(sourceOffset)%state(1,constituent)*phi
dLocalphiDot_dPhi = -sourceState(phase)%p(sourceOffset)%state(1,constituent)
end subroutine source_damage_anisobrittle_getRateAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief return array of local damage results
!--------------------------------------------------------------------------------------------------
function source_damage_anisoBrittle_postResults(phase, constituent)
use material, only: &
sourceState
2019-05-15 02:14:38 +05:30
integer, intent(in) :: &
phase, &
constituent
2019-02-13 13:35:07 +05:30
real(pReal), dimension(sum(source_damage_anisoBrittle_sizePostResult(:, &
source_damage_anisoBrittle_instance(phase)))) :: &
source_damage_anisoBrittle_postResults
2019-05-15 02:14:38 +05:30
integer :: &
instance, sourceOffset, o, c
instance = source_damage_anisoBrittle_instance(phase)
sourceOffset = source_damage_anisoBrittle_offset(phase)
2019-05-15 02:14:38 +05:30
c = 0
2019-05-15 02:14:38 +05:30
do o = 1,size(param(instance)%outputID)
2019-02-13 13:35:07 +05:30
select case(param(instance)%outputID(o))
case (damage_drivingforce_ID)
2019-05-15 02:14:38 +05:30
source_damage_anisoBrittle_postResults(c+1) = &
sourceState(phase)%p(sourceOffset)%state(1,constituent)
2019-05-15 02:14:38 +05:30
c = c + 1
end select
enddo
end function source_damage_anisoBrittle_postResults
end module source_damage_anisoBrittle