DAMASK_EICMD/src/source_damage_isoDuctile.f90

207 lines
7.9 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Luv Sharma, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine incoprorating isotropic ductile damage source mechanism
!> @details to be done
!--------------------------------------------------------------------------------------------------
module source_damage_isoDuctile
2020-02-26 23:18:33 +05:30
use prec
use debug
use IO
use discretization
use material
use config
use results
implicit none
private
integer, dimension(:), allocatable :: &
source_damage_isoDuctile_offset, & !< which source is my current damage mechanism?
source_damage_isoDuctile_instance !< instance of damage source mechanism
enum, bind(c)
enumerator :: undefined_ID, &
damage_drivingforce_ID
end enum
type, private :: tParameters !< container type for internal constitutive parameters
real(pReal) :: &
critPlasticStrain, &
N, &
aTol
integer(kind(undefined_ID)), allocatable, dimension(:) :: &
outputID
end type tParameters
type(tParameters), dimension(:), allocatable, private :: param !< containers of constitutive parameters (len Ninstance)
public :: &
source_damage_isoDuctile_init, &
source_damage_isoDuctile_dotState, &
source_damage_isoDuctile_getRateAndItsTangent, &
source_damage_isoDuctile_Results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2019-02-13 13:35:07 +05:30
subroutine source_damage_isoDuctile_init
2020-02-26 23:18:33 +05:30
integer :: Ninstance,instance,source,sourceOffset,NofMyPhase,p,i
integer(kind(undefined_ID)) :: &
outputID
character(len=pStringLen) :: &
extmsg = ''
character(len=pStringLen), dimension(:), allocatable :: &
outputs
2020-02-26 23:18:33 +05:30
write(6,'(/,a)') ' <<<+- source_'//SOURCE_DAMAGE_ISODUCTILE_LABEL//' init -+>>>'; flush(6)
2020-02-26 23:18:33 +05:30
Ninstance = count(phase_source == SOURCE_DAMAGE_ISODUCTILE_ID)
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2020-02-26 23:18:33 +05:30
allocate(source_damage_isoDuctile_offset (size(config_phase)), source=0)
allocate(source_damage_isoDuctile_instance(size(config_phase)), source=0)
allocate(param(Ninstance))
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
do p = 1, size(config_phase)
source_damage_isoDuctile_instance(p) = count(phase_source(:,1:p) == SOURCE_DAMAGE_ISODUCTILE_ID)
do source = 1, phase_Nsources(p)
if (phase_source(source,p) == SOURCE_DAMAGE_ISODUCTILE_ID) &
source_damage_isoDuctile_offset(p) = source
enddo
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
if (all(phase_source(:,p) /= SOURCE_DAMAGE_ISODUCTILE_ID)) cycle
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
associate(prm => param(source_damage_isoDuctile_instance(p)), &
config => config_phase(p))
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
prm%aTol = config%getFloat('isoductile_atol',defaultVal = 1.0e-3_pReal)
2020-02-26 23:18:33 +05:30
prm%N = config%getFloat('isoductile_ratesensitivity')
prm%critPlasticStrain = config%getFloat('isoductile_criticalplasticstrain')
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
! sanity checks
if (prm%aTol < 0.0_pReal) extmsg = trim(extmsg)//' isoductile_atol'
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
if (prm%N <= 0.0_pReal) extmsg = trim(extmsg)//' isoductile_ratesensitivity'
if (prm%critPlasticStrain <= 0.0_pReal) extmsg = trim(extmsg)//' isoductile_criticalplasticstrain'
2020-02-26 23:07:17 +05:30
2019-02-13 12:36:22 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
2020-02-26 23:18:33 +05:30
if (extmsg /= '') &
call IO_error(211,ext_msg=trim(extmsg)//'('//SOURCE_DAMAGE_ISODUCTILE_LABEL//')')
2019-02-13 12:36:22 +05:30
!--------------------------------------------------------------------------------------------------
! output pararameters
2020-02-26 23:18:33 +05:30
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
allocate(prm%outputID(0))
do i=1, size(outputs)
outputID = undefined_ID
select case(outputs(i))
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
case ('isoductile_drivingforce')
prm%outputID = [prm%outputID, damage_drivingforce_ID]
2019-02-13 13:35:07 +05:30
2020-02-26 23:18:33 +05:30
end select
2019-02-13 12:36:22 +05:30
2020-02-26 23:18:33 +05:30
enddo
2020-02-26 23:18:33 +05:30
end associate
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
NofMyPhase=count(material_phaseAt==p) * discretization_nIP
instance = source_damage_isoDuctile_instance(p)
sourceOffset = source_damage_isoDuctile_offset(p)
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
call material_allocateSourceState(p,sourceOffset,NofMyPhase,1,1,0)
sourceState(p)%p(sourceOffset)%aTolState=param(instance)%aTol
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
enddo
2020-02-26 23:07:17 +05:30
end subroutine source_damage_isoDuctile_init
2020-02-26 23:18:33 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine source_damage_isoDuctile_dotState(ipc, ip, el)
2020-02-26 23:18:33 +05:30
integer, intent(in) :: &
ipc, & !< component-ID of integration point
ip, & !< integration point
el !< element
integer :: &
phase, constituent, instance, homog, sourceOffset, damageOffset
phase = material_phaseAt(ipc,el)
constituent = material_phasememberAt(ipc,ip,el)
instance = source_damage_isoDuctile_instance(phase)
sourceOffset = source_damage_isoDuctile_offset(phase)
homog = material_homogenizationAt(el)
damageOffset = damageMapping(homog)%p(ip,el)
sourceState(phase)%p(sourceOffset)%dotState(1,constituent) = &
sum(plasticState(phase)%slipRate(:,constituent))/ &
((damage(homog)%p(damageOffset))**param(instance)%N)/ &
param(instance)%critPlasticStrain
end subroutine source_damage_isoDuctile_dotState
2020-02-26 23:07:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns local part of nonlocal damage driving force
!--------------------------------------------------------------------------------------------------
subroutine source_damage_isoDuctile_getRateAndItsTangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
2020-02-26 23:18:33 +05:30
integer, intent(in) :: &
phase, &
constituent
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
localphiDot, &
dLocalphiDot_dPhi
integer :: &
sourceOffset
2020-02-26 23:18:33 +05:30
sourceOffset = source_damage_isoDuctile_offset(phase)
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
localphiDot = 1.0_pReal &
- sourceState(phase)%p(sourceOffset)%state(1,constituent) * phi
2020-02-26 23:07:17 +05:30
2020-02-26 23:18:33 +05:30
dLocalphiDot_dPhi = -sourceState(phase)%p(sourceOffset)%state(1,constituent)
2020-02-26 23:07:17 +05:30
end subroutine source_damage_isoDuctile_getRateAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine source_damage_isoDuctile_results(phase,group)
integer, intent(in) :: phase
2020-02-26 23:07:17 +05:30
character(len=*), intent(in) :: group
integer :: sourceOffset, o, instance
2020-02-26 23:07:17 +05:30
instance = source_damage_isoDuctile_instance(phase)
sourceOffset = source_damage_isoDuctile_offset(phase)
associate(prm => param(instance), stt => sourceState(phase)%p(sourceOffset)%state)
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
case (damage_drivingforce_ID)
call results_writeDataset(group,stt,'tbd','driving force','tbd')
end select
enddo outputsLoop
end associate
end subroutine source_damage_isoDuctile_results
end module source_damage_isoDuctile