DAMASK_EICMD/code/damage_brittle.f90

431 lines
19 KiB
Fortran
Raw Normal View History

2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
! $Id$
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine incoprorating dislocation and twinning physics
!> @details to be done
!--------------------------------------------------------------------------------------------------
module damage_brittle
2014-08-01 20:24:57 +05:30
use prec, only: &
pReal, &
pInt
implicit none
private
integer(pInt), dimension(:), allocatable, public, protected :: &
damage_brittle_sizePostResults !< cumulative size of post results
2014-08-01 20:24:57 +05:30
integer(pInt), dimension(:,:), allocatable, target, public :: &
damage_brittle_sizePostResult !< size of each post result output
2014-08-01 20:24:57 +05:30
character(len=64), dimension(:,:), allocatable, target, public :: &
damage_brittle_output !< name of each post result output
2014-08-01 20:24:57 +05:30
integer(pInt), dimension(:), allocatable, target, public :: &
damage_brittle_Noutput !< number of outputs per instance of this damage
2014-08-01 20:24:57 +05:30
real(pReal), dimension(:), allocatable, private :: &
damage_brittle_aTol, &
damage_brittle_critStrainEnergy
2014-08-01 20:24:57 +05:30
enum, bind(c)
enumerator :: undefined_ID, &
local_damage_ID
end enum !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!11 ToDo
2014-08-01 20:24:57 +05:30
integer(kind(undefined_ID)), dimension(:,:), allocatable, private :: &
damage_brittle_outputID !< ID of each post result output
2014-08-01 20:24:57 +05:30
public :: &
damage_brittle_init, &
damage_brittle_stateInit, &
damage_brittle_aTolState, &
damage_brittle_dotState, &
damage_brittle_microstructure, &
constitutive_brittle_getDamage, &
constitutive_brittle_putDamage, &
damage_brittle_getDamageDiffusion33, &
damage_brittle_postResults
2014-08-01 20:24:57 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine damage_brittle_init(fileUnit)
2014-08-01 20:24:57 +05:30
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
use debug, only: &
debug_level,&
debug_constitutive,&
debug_levelBasic
use mesh, only: &
mesh_maxNips, &
mesh_NcpElems
use IO, only: &
IO_read, &
IO_lc, &
IO_getTag, &
IO_isBlank, &
IO_stringPos, &
IO_stringValue, &
IO_floatValue, &
IO_intValue, &
IO_warning, &
IO_error, &
IO_timeStamp, &
IO_EOF
use material, only: &
homogenization_maxNgrains, &
phase_damage, &
phase_damageInstance, &
phase_Noutput, &
LOCAL_DAMAGE_BRITTLE_label, &
LOCAL_DAMAGE_BRITTLE_ID, &
2014-08-01 20:24:57 +05:30
material_phase, &
damageState, &
MATERIAL_partPhase
use numerics,only: &
worldrank, &
2014-08-01 20:24:57 +05:30
numerics_integrator
implicit none
integer(pInt), intent(in) :: fileUnit
integer(pInt), parameter :: MAXNCHUNKS = 7_pInt
integer(pInt), dimension(1+2*MAXNCHUNKS) :: positions
integer(pInt) :: maxNinstance,mySize=0_pInt,phase,instance,o
integer(pInt) :: sizeState, sizeDotState
integer(pInt) :: NofMyPhase
character(len=65536) :: &
tag = '', &
line = ''
mainProcess: if (worldrank == 0) then
write(6,'(/,a)') ' <<<+- damage_'//LOCAL_DAMAGE_BRITTLE_label//' init -+>>>'
write(6,'(a)') ' $Id$'
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
2014-08-01 20:24:57 +05:30
#include "compilation_info.f90"
endif mainProcess
maxNinstance = int(count(phase_damage == LOCAL_DAMAGE_BRITTLE_ID),pInt)
2014-08-01 20:24:57 +05:30
if (maxNinstance == 0_pInt) return
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0_pInt) &
write(6,'(a16,1x,i5,/)') '# instances:',maxNinstance
allocate(damage_brittle_sizePostResults(maxNinstance), source=0_pInt)
allocate(damage_brittle_sizePostResult(maxval(phase_Noutput),maxNinstance),source=0_pInt)
allocate(damage_brittle_output(maxval(phase_Noutput),maxNinstance))
damage_brittle_output = ''
allocate(damage_brittle_outputID(maxval(phase_Noutput),maxNinstance), source=undefined_ID)
allocate(damage_brittle_Noutput(maxNinstance), source=0_pInt)
allocate(damage_brittle_critStrainEnergy(maxNinstance), source=0.0_pReal)
allocate(damage_brittle_aTol(maxNinstance), source=0.0_pReal)
2014-08-01 20:24:57 +05:30
rewind(fileUnit)
phase = 0_pInt
do while (trim(line) /= IO_EOF .and. IO_lc(IO_getTag(line,'<','>')) /= MATERIAL_partPhase) ! wind forward to <phase>
line = IO_read(fileUnit)
enddo
parsingFile: do while (trim(line) /= IO_EOF) ! read through sections of phase part
line = IO_read(fileUnit)
if (IO_isBlank(line)) cycle ! skip empty lines
if (IO_getTag(line,'<','>') /= '') then ! stop at next part
line = IO_read(fileUnit, .true.) ! reset IO_read
exit
endif
if (IO_getTag(line,'[',']') /= '') then ! next phase section
phase = phase + 1_pInt ! advance phase section counter
cycle ! skip to next line
endif
if (phase > 0_pInt ) then; if (phase_damage(phase) == LOCAL_DAMAGE_BRITTLE_ID) then ! do not short-circuit here (.and. with next if statemen). It's not safe in Fortran
2014-08-01 20:24:57 +05:30
instance = phase_damageInstance(phase) ! which instance of my damage is present phase
positions = IO_stringPos(line,MAXNCHUNKS)
tag = IO_lc(IO_stringValue(line,positions,1_pInt)) ! extract key
select case(tag)
case ('(output)')
select case(IO_lc(IO_stringValue(line,positions,2_pInt)))
case ('local_damage')
damage_brittle_Noutput(instance) = damage_brittle_Noutput(instance) + 1_pInt
damage_brittle_outputID(damage_brittle_Noutput(instance),instance) = local_damage_ID
damage_brittle_output(damage_brittle_Noutput(instance),instance) = &
2014-08-01 20:24:57 +05:30
IO_lc(IO_stringValue(line,positions,2_pInt))
end select
case ('critical_strain_energy')
damage_brittle_critStrainEnergy(instance) = IO_floatValue(line,positions,2_pInt)
2014-08-01 20:24:57 +05:30
case ('atol_damage')
damage_brittle_aTol(instance) = IO_floatValue(line,positions,2_pInt)
2014-08-01 20:24:57 +05:30
end select
endif; endif
enddo parsingFile
initializeInstances: do phase = 1_pInt, size(phase_damage)
if (phase_damage(phase) == LOCAL_DAMAGE_BRITTLE_ID) then
2014-08-01 20:24:57 +05:30
NofMyPhase=count(material_phase==phase)
instance = phase_damageInstance(phase)
!--------------------------------------------------------------------------------------------------
! Determine size of postResults array
outputsLoop: do o = 1_pInt,damage_brittle_Noutput(instance)
select case(damage_brittle_outputID(o,instance))
2014-08-01 20:24:57 +05:30
case(local_damage_ID)
mySize = 1_pInt
end select
if (mySize > 0_pInt) then ! any meaningful output found
damage_brittle_sizePostResult(o,instance) = mySize
damage_brittle_sizePostResults(instance) = damage_brittle_sizePostResults(instance) + mySize
2014-08-01 20:24:57 +05:30
endif
enddo outputsLoop
! Determine size of state array
sizeDotState = 1_pInt
2014-08-01 20:24:57 +05:30
sizeState = 2_pInt
damageState(phase)%sizeState = sizeState
damageState(phase)%sizeDotState = sizeDotState
damageState(phase)%sizePostResults = damage_brittle_sizePostResults(instance)
2014-08-01 20:24:57 +05:30
allocate(damageState(phase)%aTolState (sizeState), source=0.0_pReal)
allocate(damageState(phase)%state0 (sizeState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%partionedState0 (sizeState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%subState0 (sizeState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%state (sizeState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%state_backup (sizeState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%dotState (sizeDotState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%deltaState (sizeDotState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%dotState_backup (sizeDotState,NofMyPhase), source=0.0_pReal)
if (any(numerics_integrator == 1_pInt)) then
allocate(damageState(phase)%previousDotState (sizeDotState,NofMyPhase), source=0.0_pReal)
allocate(damageState(phase)%previousDotState2 (sizeDotState,NofMyPhase), source=0.0_pReal)
endif
if (any(numerics_integrator == 4_pInt)) &
allocate(damageState(phase)%RK4dotState (sizeDotState,NofMyPhase), source=0.0_pReal)
if (any(numerics_integrator == 5_pInt)) &
allocate(damageState(phase)%RKCK45dotState (6,sizeDotState,NofMyPhase),source=0.0_pReal)
2014-10-13 23:22:33 +05:30
call damage_brittle_stateInit(phase)
call damage_brittle_aTolState(phase,instance)
2014-08-01 20:24:57 +05:30
endif
enddo initializeInstances
end subroutine damage_brittle_init
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief sets the relevant NEW state values for a given instance of this damage
!--------------------------------------------------------------------------------------------------
2014-10-13 23:22:33 +05:30
subroutine damage_brittle_stateInit(phase)
2014-08-01 20:24:57 +05:30
use material, only: &
damageState
implicit none
integer(pInt), intent(in) :: phase !< number specifying the phase of the damage
real(pReal), dimension(damageState(phase)%sizeState) :: tempState
tempState(1) = 1.0_pReal
tempState(2) = 1.0_pReal
2014-08-01 20:24:57 +05:30
damageState(phase)%state = spread(tempState,2,size(damageState(phase)%state(1,:)))
damageState(phase)%state0 = damageState(phase)%state
damageState(phase)%partionedState0 = damageState(phase)%state
end subroutine damage_brittle_stateInit
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief sets the relevant state values for a given instance of this damage
!--------------------------------------------------------------------------------------------------
subroutine damage_brittle_aTolState(phase,instance)
2014-08-01 20:24:57 +05:30
use material, only: &
damageState
implicit none
integer(pInt), intent(in) :: &
phase, &
instance ! number specifying the current instance of the damage
real(pReal), dimension(damageState(phase)%sizeState) :: tempTol
tempTol = damage_brittle_aTol(instance)
2014-08-01 20:24:57 +05:30
damageState(phase)%aTolState = tempTol
end subroutine damage_brittle_aTolState
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine damage_brittle_dotState(ipc, ip, el)
use material, only: &
mappingConstitutive, &
phase_damageInstance, &
damageState
use lattice, only: &
lattice_DamageMobility
implicit none
integer(pInt), intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
integer(pInt) :: &
phase, constituent
phase = mappingConstitutive(2,ipc,ip,el)
constituent = mappingConstitutive(1,ipc,ip,el)
damageState(phase)%dotState(1,constituent) = &
(1.0_pReal/lattice_DamageMobility(phase))* &
(damageState(phase)%state(2,constituent) - &
damageState(phase)%state(1,constituent))
end subroutine damage_brittle_dotState
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine damage_brittle_microstructure(Tstar_v, Fe, ipc, ip, el)
2014-08-01 20:24:57 +05:30
use material, only: &
mappingConstitutive, &
phase_damageInstance, &
2014-08-01 20:24:57 +05:30
damageState
use math, only: &
math_Mandel6to33, &
2014-08-01 20:24:57 +05:30
math_mul33x33, &
math_transpose33, &
math_I3
implicit none
integer(pInt), intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
real(pReal), intent(in), dimension(6) :: &
Tstar_v !< 2nd Piola Kirchhoff stress tensor (Mandel)
real(pReal), intent(in), dimension(3,3) :: &
Fe
integer(pInt) :: &
phase, constituent, instance
2014-08-01 20:24:57 +05:30
real(pReal) :: &
2014-10-06 22:31:39 +05:30
strain(3,3)
2014-08-01 20:24:57 +05:30
phase = mappingConstitutive(2,ipc,ip,el)
constituent = mappingConstitutive(1,ipc,ip,el)
instance = phase_damageInstance(phase)
strain = 0.5_pReal*(math_mul33x33(math_transpose33(Fe),Fe)-math_I3)
damageState(phase)%state(2,constituent) = min(damageState(phase)%state0(2,constituent), &
damage_brittle_critStrainEnergy(instance)/ &
sum(abs(math_Mandel6to33(Tstar_v)*strain)))
end subroutine damage_brittle_microstructure
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns temperature based on local damage model state layout
!--------------------------------------------------------------------------------------------------
function constitutive_brittle_getDamage(ipc, ip, el)
use material, only: &
mappingConstitutive, &
damageState
implicit none
integer(pInt), intent(in) :: &
ipc, & !< grain number
ip, & !< integration point number
el !< element number
real(pReal) :: constitutive_brittle_getDamage
2014-08-01 20:24:57 +05:30
constitutive_brittle_getDamage = &
damageState(mappingConstitutive(2,ipc,ip,el))%state(1,mappingConstitutive(1,ipc,ip,el))
end function constitutive_brittle_getDamage
!--------------------------------------------------------------------------------------------------
!> @brief returns temperature based on local damage model state layout
!--------------------------------------------------------------------------------------------------
subroutine constitutive_brittle_putDamage(ipc, ip, el, localDamage)
use material, only: &
mappingConstitutive, &
damageState
implicit none
integer(pInt), intent(in) :: &
ipc, & !< grain number
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: localDamage
damageState(mappingConstitutive(2,ipc,ip,el))%state(1,mappingConstitutive(1,ipc,ip,el)) = &
localDamage
end subroutine constitutive_brittle_putDamage
!--------------------------------------------------------------------------------------------------
!> @brief returns brittle damage diffusion tensor
!--------------------------------------------------------------------------------------------------
function damage_brittle_getDamageDiffusion33(ipc, ip, el)
use lattice, only: &
lattice_DamageDiffusion33
use material, only: &
mappingConstitutive, &
damageState
implicit none
integer(pInt), intent(in) :: &
ipc, & !< grain number
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
damage_brittle_getDamageDiffusion33
integer(pInt) :: &
phase, constituent
phase = mappingConstitutive(2,ipc,ip,el)
constituent = mappingConstitutive(1,ipc,ip,el)
damage_brittle_getDamageDiffusion33 = &
damageState(phase)%state(2,constituent)* &
lattice_DamageDiffusion33(1:3,1:3,phase)
end function damage_brittle_getDamageDiffusion33
2014-08-01 20:24:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results
!--------------------------------------------------------------------------------------------------
function damage_brittle_postResults(ipc,ip,el)
2014-08-01 20:24:57 +05:30
use material, only: &
mappingConstitutive, &
phase_damageInstance,&
damageState
implicit none
integer(pInt), intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
real(pReal), dimension(damage_brittle_sizePostResults(phase_damageInstance(mappingConstitutive(2,ipc,ip,el)))) :: &
damage_brittle_postResults
2014-08-01 20:24:57 +05:30
integer(pInt) :: &
instance, phase, constituent, o, c
phase = mappingConstitutive(2,ipc,ip,el)
constituent = mappingConstitutive(1,ipc,ip,el)
instance = phase_damageInstance(phase)
c = 0_pInt
damage_brittle_postResults = 0.0_pReal
2014-08-01 20:24:57 +05:30
do o = 1_pInt,damage_brittle_Noutput(instance)
select case(damage_brittle_outputID(o,instance))
2014-08-01 20:24:57 +05:30
case (local_damage_ID)
damage_brittle_postResults(c+1_pInt) = damageState(phase)%state(1,constituent)
2014-08-01 20:24:57 +05:30
c = c + 1
end select
enddo
end function damage_brittle_postResults
2014-08-01 20:24:57 +05:30
end module damage_brittle