2015-05-28 22:32:23 +05:30
!--------------------------------------------------------------------------------------------------
!> @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 ductile damage source mechanism
!> @details to be done
!--------------------------------------------------------------------------------------------------
module source_damage_anisoDuctile
use prec , only : &
pReal , &
pInt
implicit none
private
integer ( pInt ) , dimension ( : ) , allocatable , public , protected :: &
source_damage_anisoDuctile_sizePostResults , & !< cumulative size of post results
source_damage_anisoDuctile_offset , & !< which source is my current damage mechanism?
source_damage_anisoDuctile_instance !< instance of damage source mechanism
integer ( pInt ) , dimension ( : , : ) , allocatable , target , public :: &
source_damage_anisoDuctile_sizePostResult !< size of each post result output
character ( len = 64 ) , dimension ( : , : ) , allocatable , target , public :: &
source_damage_anisoDuctile_output !< name of each post result output
integer ( pInt ) , dimension ( : ) , allocatable , target , public :: &
source_damage_anisoDuctile_Noutput !< number of outputs per instance of this damage
integer ( pInt ) , dimension ( : ) , allocatable , private :: &
source_damage_anisoDuctile_totalNslip !< total number of slip systems
integer ( pInt ) , dimension ( : , : ) , allocatable , private :: &
source_damage_anisoDuctile_Nslip !< number of slip systems per family
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
source_damage_anisoDuctile_critPlasticStrain
real ( pReal ) , dimension ( : , : ) , allocatable , private :: &
source_damage_anisoDuctile_critLoad
enum , bind ( c )
enumerator :: undefined_ID , &
damage_drivingforce_ID
end enum
integer ( kind ( undefined_ID ) ) , dimension ( : , : ) , allocatable , private :: &
source_damage_anisoDuctile_outputID !< ID of each post result output
2018-12-31 01:38:48 +05:30
type , private :: tParameters !< container type for internal constitutive parameters
real ( pReal ) :: &
aTol , &
N
real ( pReal ) , dimension ( : ) , allocatable :: &
critPlasticStrain , &
critLoad
integer ( pInt ) :: &
totalNslip
integer ( pInt ) , dimension ( : ) , allocatable :: &
Nslip
2019-02-13 12:36:22 +05:30
integer ( kind ( undefined_ID ) ) , allocatable , dimension ( : ) :: &
outputID
2018-12-31 01:38:48 +05:30
end type tParameters
type ( tParameters ) , dimension ( : ) , allocatable , private :: param !< containers of constitutive parameters (len Ninstance)
2015-05-28 22:32:23 +05:30
public :: &
source_damage_anisoDuctile_init , &
source_damage_anisoDuctile_dotState , &
source_damage_anisoDuctile_getRateAndItsTangent , &
source_damage_anisoDuctile_postResults
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisoDuctile_init ( fileUnit )
2018-02-02 17:06:09 +05:30
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
2017-10-05 20:05:34 +05:30
use , intrinsic :: iso_fortran_env , only : &
compiler_version , &
compiler_options
#endif
2019-02-13 10:58:29 +05:30
use prec , only : &
pStringLen
2015-05-28 22:32:23 +05:30
use debug , only : &
debug_level , &
debug_constitutive , &
debug_levelBasic
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 : &
2019-02-13 11:52:37 +05:30
material_allocateSourceState , &
2015-05-28 22:32:23 +05:30
phase_source , &
phase_Nsources , &
phase_Noutput , &
SOURCE_damage_anisoDuctile_label , &
SOURCE_damage_anisoDuctile_ID , &
material_phase , &
2018-06-10 21:31:52 +05:30
sourceState
2018-06-14 10:09:49 +05:30
use config , only : &
2018-12-31 01:38:48 +05:30
config_phase , &
2018-06-10 21:31:52 +05:30
material_Nphase , &
2015-05-28 22:32:23 +05:30
MATERIAL_partPhase
use lattice , only : &
lattice_maxNslipFamily , &
lattice_NslipSystem
implicit none
integer ( pInt ) , intent ( in ) :: fileUnit
2015-08-28 13:08:48 +05:30
integer ( pInt ) , allocatable , dimension ( : ) :: chunkPos
2018-12-31 01:38:48 +05:30
integer ( pInt ) :: Ninstance , mySize = 0_pInt , phase , instance , source , sourceOffset , o
2019-02-13 12:36:22 +05:30
integer ( pInt ) :: NofMyPhase , p , i
integer ( pInt ) :: Nchunks_SlipFamilies = 0_pInt , j
2015-05-28 22:32:23 +05:30
character ( len = 65536 ) :: &
tag = '' , &
line = ''
2019-02-13 10:58:29 +05:30
integer ( pInt ) , dimension ( 0 ) , parameter :: emptyIntArray = [ integer ( pInt ) :: ]
2019-02-13 12:36:22 +05:30
character ( len = 65536 ) , dimension ( 0 ) , parameter :: emptyStringArray = [ character ( len = 65536 ) :: ]
integer ( kind ( undefined_ID ) ) :: &
outputID
character ( len = pStringLen ) :: &
extmsg = ''
character ( len = 65536 ) , dimension ( : ) , allocatable :: &
outputs
2015-05-28 22:32:23 +05:30
2019-02-13 12:36:22 +05:30
write ( 6 , '(/,a)' ) ' <<<+- source_' / / SOURCE_DAMAGE_ANISODUCTILE_LABEL / / ' init -+>>>'
2016-07-25 23:42:00 +05:30
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2015-05-28 22:32:23 +05:30
#include "compilation_info.f90"
2018-12-31 01:38:48 +05:30
Ninstance = int ( count ( phase_source == SOURCE_damage_anisoDuctile_ID ) , pInt )
if ( Ninstance == 0_pInt ) return
2015-05-28 22:32:23 +05:30
if ( iand ( debug_level ( debug_constitutive ) , debug_levelBasic ) / = 0_pInt ) &
2018-12-31 01:38:48 +05:30
write ( 6 , '(a16,1x,i5,/)' ) '# instances:' , Ninstance
2015-05-28 22:32:23 +05:30
allocate ( source_damage_anisoDuctile_offset ( material_Nphase ) , source = 0_pInt )
allocate ( source_damage_anisoDuctile_instance ( material_Nphase ) , source = 0_pInt )
do phase = 1 , material_Nphase
source_damage_anisoDuctile_instance ( phase ) = count ( phase_source ( : , 1 : phase ) == source_damage_anisoDuctile_ID )
do source = 1 , phase_Nsources ( phase )
if ( phase_source ( source , phase ) == source_damage_anisoDuctile_ID ) &
source_damage_anisoDuctile_offset ( phase ) = source
enddo
enddo
2018-12-31 01:38:48 +05:30
allocate ( source_damage_anisoDuctile_sizePostResults ( Ninstance ) , source = 0_pInt )
allocate ( source_damage_anisoDuctile_sizePostResult ( maxval ( phase_Noutput ) , Ninstance ) , source = 0_pInt )
allocate ( source_damage_anisoDuctile_output ( maxval ( phase_Noutput ) , Ninstance ) )
2015-05-28 22:32:23 +05:30
source_damage_anisoDuctile_output = ''
2018-12-31 01:38:48 +05:30
allocate ( source_damage_anisoDuctile_outputID ( maxval ( phase_Noutput ) , Ninstance ) , source = undefined_ID )
allocate ( source_damage_anisoDuctile_Noutput ( Ninstance ) , source = 0_pInt )
allocate ( source_damage_anisoDuctile_critLoad ( lattice_maxNslipFamily , Ninstance ) , source = 0.0_pReal )
allocate ( source_damage_anisoDuctile_critPlasticStrain ( lattice_maxNslipFamily , Ninstance ) , source = 0.0_pReal )
allocate ( source_damage_anisoDuctile_Nslip ( lattice_maxNslipFamily , Ninstance ) , source = 0_pInt )
2019-02-13 10:58:29 +05:30
allocate ( source_damage_anisoDuctile_totalNslip ( Ninstance ) , source = 0_pInt )
2018-12-31 01:38:48 +05:30
2019-02-13 10:58:29 +05:30
allocate ( param ( Ninstance ) )
2018-12-31 01:38:48 +05:30
do p = 1 , size ( config_phase )
2019-02-13 10:58:29 +05:30
if ( all ( phase_source ( : , p ) / = SOURCE_DAMAGE_ANISODUCTILE_ID ) ) cycle
associate ( prm = > param ( source_damage_anisoDuctile_instance ( p ) ) , &
config = > config_phase ( p ) )
prm % aTol = config % getFloat ( 'anisoductile_atol' , defaultVal = 1.0e-3_pReal )
prm % N = config % getFloat ( 'anisoductile_ratesensitivity' )
! sanity checks
if ( prm % aTol < 0.0_pReal ) extmsg = trim ( extmsg ) / / ' anisoductile_atol'
if ( prm % N < = 0.0_pReal ) extmsg = trim ( extmsg ) / / ' anisoductile_ratesensitivity'
prm % Nslip = config % getInts ( 'nslip' , defaultVal = emptyIntArray )
2019-02-13 12:36:22 +05:30
!--------------------------------------------------------------------------------------------------
! exit if any parameter is out of range
if ( extmsg / = '' ) &
call IO_error ( 211_pInt , ext_msg = trim ( extmsg ) / / '(' / / SOURCE_DAMAGE_ANISODUCTILE_LABEL / / ')' )
!--------------------------------------------------------------------------------------------------
! output pararameters
outputs = config % getStrings ( '(output)' , defaultVal = emptyStringArray )
allocate ( prm % outputID ( 0 ) )
do i = 1_pInt , size ( outputs )
outputID = undefined_ID
select case ( outputs ( i ) )
case ( 'anisoductile_drivingforce' )
end select
enddo
2019-02-13 10:58:29 +05:30
end associate
2018-12-31 01:38:48 +05:30
enddo
2015-05-28 22:32:23 +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 ( any ( phase_source ( : , phase ) == SOURCE_damage_anisoDuctile_ID ) ) then ! do not short-circuit here (.and. with next if statemen). It's not safe in Fortran
instance = source_damage_anisoDuctile_instance ( phase ) ! which instance of my damage is present phase
2015-08-28 13:08:48 +05:30
chunkPos = IO_stringPos ( line )
tag = IO_lc ( IO_stringValue ( line , chunkPos , 1_pInt ) ) ! extract key
2015-05-28 22:32:23 +05:30
select case ( tag )
case ( '(output)' )
2015-08-28 13:08:48 +05:30
select case ( IO_lc ( IO_stringValue ( line , chunkPos , 2_pInt ) ) )
2015-05-28 22:32:23 +05:30
case ( 'anisoductile_drivingforce' )
source_damage_anisoDuctile_Noutput ( instance ) = source_damage_anisoDuctile_Noutput ( instance ) + 1_pInt
source_damage_anisoDuctile_outputID ( source_damage_anisoDuctile_Noutput ( instance ) , instance ) = damage_drivingforce_ID
source_damage_anisoDuctile_output ( source_damage_anisoDuctile_Noutput ( instance ) , instance ) = &
2015-08-28 13:08:48 +05:30
IO_lc ( IO_stringValue ( line , chunkPos , 2_pInt ) )
2015-05-28 22:32:23 +05:30
end select
case ( 'nslip' ) !
2015-08-28 13:08:48 +05:30
Nchunks_SlipFamilies = chunkPos ( 1 ) - 1_pInt
2015-05-28 22:32:23 +05:30
do j = 1_pInt , Nchunks_SlipFamilies
2015-08-28 13:08:48 +05:30
source_damage_anisoDuctile_Nslip ( j , instance ) = IO_intValue ( line , chunkPos , 1_pInt + j )
2015-05-28 22:32:23 +05:30
enddo
case ( 'anisoductile_criticalplasticstrain' )
do j = 1_pInt , Nchunks_SlipFamilies
2015-08-28 13:08:48 +05:30
source_damage_anisoDuctile_critPlasticStrain ( j , instance ) = IO_floatValue ( line , chunkPos , 1_pInt + j )
2015-05-28 22:32:23 +05:30
enddo
case ( 'anisoductile_criticalload' )
do j = 1_pInt , Nchunks_SlipFamilies
2015-08-28 13:08:48 +05:30
source_damage_anisoDuctile_critLoad ( j , instance ) = IO_floatValue ( line , chunkPos , 1_pInt + j )
2015-05-28 22:32:23 +05:30
enddo
end select
endif ; endif
enddo parsingFile
!--------------------------------------------------------------------------------------------------
! sanity checks
sanityChecks : do phase = 1_pInt , size ( phase_source )
myPhase : if ( any ( phase_source ( : , phase ) == SOURCE_damage_anisoDuctile_ID ) ) then
instance = source_damage_anisoDuctile_instance ( phase )
source_damage_anisoDuctile_Nslip ( 1 : lattice_maxNslipFamily , instance ) = &
min ( lattice_NslipSystem ( 1 : lattice_maxNslipFamily , phase ) , & ! limit active cleavage systems per family to min of available and requested
source_damage_anisoDuctile_Nslip ( 1 : lattice_maxNslipFamily , instance ) )
source_damage_anisoDuctile_totalNslip ( instance ) = sum ( source_damage_anisoDuctile_Nslip ( : , instance ) )
2019-02-13 10:58:29 +05:30
2015-05-28 22:32:23 +05:30
if ( any ( source_damage_anisoDuctile_critPlasticStrain ( : , instance ) < 0.0_pReal ) ) &
call IO_error ( 211_pInt , el = instance , ext_msg = 'criticaPlasticStrain (' / / SOURCE_damage_anisoDuctile_LABEL / / ')' )
2019-02-13 10:58:29 +05:30
2015-05-28 22:32:23 +05:30
endif myPhase
enddo sanityChecks
initializeInstances : do phase = 1_pInt , material_Nphase
if ( any ( phase_source ( : , phase ) == SOURCE_damage_anisoDuctile_ID ) ) then
NofMyPhase = count ( material_phase == phase )
instance = source_damage_anisoDuctile_instance ( phase )
sourceOffset = source_damage_anisoDuctile_offset ( phase )
!--------------------------------------------------------------------------------------------------
! Determine size of postResults array
outputsLoop : do o = 1_pInt , source_damage_anisoDuctile_Noutput ( instance )
select case ( source_damage_anisoDuctile_outputID ( o , instance ) )
case ( damage_drivingforce_ID )
mySize = 1_pInt
end select
if ( mySize > 0_pInt ) then ! any meaningful output found
source_damage_anisoDuctile_sizePostResult ( o , instance ) = mySize
source_damage_anisoDuctile_sizePostResults ( instance ) = source_damage_anisoDuctile_sizePostResults ( instance ) + mySize
endif
enddo outputsLoop
2019-02-13 11:52:37 +05:30
call material_allocateSourceState ( phase , sourceOffset , NofMyPhase , 1_pInt )
2015-05-28 22:32:23 +05:30
sourceState ( phase ) % p ( sourceOffset ) % sizePostResults = source_damage_anisoDuctile_sizePostResults ( instance )
2019-02-13 11:52:37 +05:30
sourceState ( phase ) % p ( sourceOffset ) % aTolState = param ( instance ) % aTol
2015-05-28 22:32:23 +05:30
endif
enddo initializeInstances
end subroutine source_damage_anisoDuctile_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates derived quantities from state
!--------------------------------------------------------------------------------------------------
subroutine source_damage_anisoDuctile_dotState ( ipc , ip , el )
use material , only : &
2016-01-15 05:49:44 +05:30
phaseAt , phasememberAt , &
2015-05-28 22:32:23 +05:30
plasticState , &
sourceState , &
material_homog , &
damage , &
damageMapping
use lattice , only : &
lattice_maxNslipFamily
implicit none
integer ( pInt ) , intent ( in ) :: &
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
integer ( pInt ) :: &
phase , &
constituent , &
sourceOffset , &
homog , damageOffset , &
instance , &
index , f , i
2016-01-15 05:49:44 +05:30
phase = phaseAt ( ipc , ip , el )
constituent = phasememberAt ( ipc , ip , el )
2015-05-28 22:32:23 +05:30
instance = source_damage_anisoDuctile_instance ( phase )
sourceOffset = source_damage_anisoDuctile_offset ( phase )
homog = material_homog ( ip , el )
damageOffset = damageMapping ( homog ) % p ( ip , el )
index = 1_pInt
sourceState ( phase ) % p ( sourceOffset ) % dotState ( 1 , constituent ) = 0.0_pReal
do f = 1_pInt , lattice_maxNslipFamily
do i = 1_pInt , source_damage_anisoDuctile_Nslip ( f , instance ) ! process each (active) slip system in family
sourceState ( phase ) % p ( sourceOffset ) % dotState ( 1 , constituent ) = &
sourceState ( phase ) % p ( sourceOffset ) % dotState ( 1 , constituent ) + &
plasticState ( phase ) % slipRate ( index , constituent ) / &
2019-02-13 10:58:29 +05:30
( ( damage ( homog ) % p ( damageOffset ) ) ** param ( instance ) % N ) / &
2015-05-28 22:32:23 +05:30
source_damage_anisoDuctile_critPlasticStrain ( f , instance )
index = index + 1_pInt
enddo
enddo
end subroutine source_damage_anisoDuctile_dotState
!--------------------------------------------------------------------------------------------------
!> @brief returns local part of nonlocal damage driving force
!--------------------------------------------------------------------------------------------------
2018-12-31 01:38:48 +05:30
subroutine source_damage_anisoDuctile_getRateAndItsTangent ( localphiDot , dLocalphiDot_dPhi , phi , phase , constituent )
2015-05-28 22:32:23 +05:30
use material , only : &
sourceState
implicit none
integer ( pInt ) , intent ( in ) :: &
2018-12-31 01:38:48 +05:30
phase , &
constituent
2015-05-28 22:32:23 +05:30
real ( pReal ) , intent ( in ) :: &
phi
real ( pReal ) , intent ( out ) :: &
localphiDot , &
dLocalphiDot_dPhi
integer ( pInt ) :: &
2018-12-31 01:38:48 +05:30
sourceOffset
2015-05-28 22:32:23 +05:30
sourceOffset = source_damage_anisoDuctile_offset ( phase )
localphiDot = 1.0_pReal - &
2015-06-11 14:33:51 +05:30
sourceState ( phase ) % p ( sourceOffset ) % state ( 1 , constituent ) * &
2015-05-28 22:32:23 +05:30
phi
2015-06-11 14:33:51 +05:30
dLocalphiDot_dPhi = - sourceState ( phase ) % p ( sourceOffset ) % state ( 1 , constituent )
2015-05-28 22:32:23 +05:30
end subroutine source_damage_anisoDuctile_getRateAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief return array of local damage results
!--------------------------------------------------------------------------------------------------
2018-12-31 01:38:48 +05:30
function source_damage_anisoDuctile_postResults ( phase , constituent )
2015-05-28 22:32:23 +05:30
use material , only : &
sourceState
implicit none
2018-12-31 01:38:48 +05:30
integer ( pInt ) , intent ( in ) :: &
phase , &
constituent
2015-05-28 22:32:23 +05:30
real ( pReal ) , dimension ( source_damage_anisoDuctile_sizePostResults ( &
2018-12-31 01:38:48 +05:30
source_damage_anisoDuctile_instance ( phase ) ) ) :: &
2015-05-28 22:32:23 +05:30
source_damage_anisoDuctile_postResults
integer ( pInt ) :: &
2018-12-31 01:38:48 +05:30
instance , sourceOffset , o , c
2015-05-28 22:32:23 +05:30
instance = source_damage_anisoDuctile_instance ( phase )
sourceOffset = source_damage_anisoDuctile_offset ( phase )
c = 0_pInt
source_damage_anisoDuctile_postResults = 0.0_pReal
do o = 1_pInt , source_damage_anisoDuctile_Noutput ( instance )
select case ( source_damage_anisoDuctile_outputID ( o , instance ) )
case ( damage_drivingforce_ID )
source_damage_anisoDuctile_postResults ( c + 1_pInt ) = &
sourceState ( phase ) % p ( sourceOffset ) % state ( 1 , constituent )
c = c + 1_pInt
end select
enddo
end function source_damage_anisoDuctile_postResults
end module source_damage_anisoDuctile