2012-10-02 18:23:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
2013-03-06 20:11:15 +05:30
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
2012-10-02 18:23:25 +05:30
!> @brief elasticity, plasticity, internal microstructure state
!--------------------------------------------------------------------------------------------------
module constitutive
2013-01-16 15:44:57 +05:30
use prec , only : &
2014-07-02 17:57:39 +05:30
pInt
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
implicit none
private
2014-05-27 20:16:03 +05:30
integer ( pInt ) , public , protected :: &
2015-05-28 22:32:23 +05:30
constitutive_plasticity_maxSizePostResults , &
constitutive_plasticity_maxSizeDotState , &
constitutive_source_maxSizePostResults , &
constitutive_source_maxSizeDotState
2014-06-26 19:23:12 +05:30
2015-10-14 00:22:01 +05:30
public :: &
2013-01-16 15:44:57 +05:30
constitutive_init , &
constitutive_homogenizedC , &
constitutive_microstructure , &
2018-08-28 18:37:39 +05:30
constitutive_LpAndItsTangents , &
constitutive_LiAndItsTangents , &
2015-07-24 20:17:18 +05:30
constitutive_initialFi , &
2018-08-28 18:24:36 +05:30
constitutive_SandItsTangents , &
2013-01-16 15:44:57 +05:30
constitutive_collectDotState , &
constitutive_collectDeltaState , &
2018-12-05 04:25:39 +05:30
constitutive_postResults , &
constitutive_results
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
private :: &
2018-08-28 18:24:36 +05:30
constitutive_hooke_SandItsTangents
2015-10-14 00:22:01 +05:30
2012-03-09 01:55:28 +05:30
contains
2012-10-02 18:23:25 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief allocates arrays pointing to array of the various constitutive modules
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
subroutine constitutive_init ( )
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
2014-07-02 17:57:39 +05:30
use prec , only : &
pReal
2013-01-16 15:44:57 +05:30
use debug , only : &
debug_constitutive , &
debug_levelBasic
use numerics , only : &
2015-04-21 21:41:30 +05:30
worldrank
2013-01-16 15:44:57 +05:30
use IO , only : &
IO_error , &
IO_open_file , &
2016-01-16 22:57:19 +05:30
IO_checkAndRewind , &
2013-01-16 15:44:57 +05:30
IO_open_jobFile_stat , &
IO_write_jobFile , &
2013-09-18 19:37:55 +05:30
IO_write_jobIntFile , &
2013-02-25 22:04:59 +05:30
IO_timeStamp
2018-06-27 00:01:52 +05:30
use config , only : &
2018-08-22 15:21:23 +05:30
config_phase
2013-01-16 15:44:57 +05:30
use mesh , only : &
FE_geomtype
2018-06-14 10:09:49 +05:30
use config , only : &
2013-01-16 15:44:57 +05:30
material_Nphase , &
2015-10-14 00:22:01 +05:30
material_localFileExt , &
2013-01-16 15:44:57 +05:30
phase_name , &
2018-08-22 18:00:51 +05:30
material_configFile , &
config_deallocate
2018-06-10 21:31:52 +05:30
use material , only : &
material_phase , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
phase_plasticityInstance , &
2015-05-28 22:32:23 +05:30
phase_Nsources , &
phase_source , &
phase_kinematics , &
2014-10-09 19:38:32 +05:30
ELASTICITY_hooke_ID , &
PLASTICITY_none_ID , &
2016-01-09 01:15:20 +05:30
PLASTICITY_isotropic_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_phenopowerlaw_ID , &
2017-09-30 04:06:28 +05:30
PLASTICITY_kinehardening_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_dislotwin_ID , &
2015-01-15 16:26:15 +05:30
PLASTICITY_disloucla_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_nonlocal_ID , &
2015-05-28 22:32:23 +05:30
SOURCE_thermal_dissipation_ID , &
2015-07-27 16:39:37 +05:30
SOURCE_thermal_externalheat_ID , &
2015-05-28 22:32:23 +05:30
SOURCE_damage_isoBrittle_ID , &
SOURCE_damage_isoDuctile_ID , &
SOURCE_damage_anisoBrittle_ID , &
SOURCE_damage_anisoDuctile_ID , &
KINEMATICS_cleavage_opening_ID , &
KINEMATICS_slipplane_opening_ID , &
KINEMATICS_thermal_expansion_ID , &
2013-11-27 13:34:05 +05:30
ELASTICITY_HOOKE_label , &
PLASTICITY_NONE_label , &
2016-01-09 01:15:20 +05:30
PLASTICITY_ISOTROPIC_label , &
2013-11-27 13:34:05 +05:30
PLASTICITY_PHENOPOWERLAW_label , &
2017-09-30 04:06:28 +05:30
PLASTICITY_KINEHARDENING_label , &
2013-11-27 13:34:05 +05:30
PLASTICITY_DISLOTWIN_label , &
2015-01-15 16:26:15 +05:30
PLASTICITY_DISLOUCLA_label , &
2014-08-08 16:34:40 +05:30
PLASTICITY_NONLOCAL_label , &
2015-05-28 22:32:23 +05:30
SOURCE_thermal_dissipation_label , &
2015-07-27 16:39:37 +05:30
SOURCE_thermal_externalheat_label , &
2015-05-28 22:32:23 +05:30
SOURCE_damage_isoBrittle_label , &
SOURCE_damage_isoDuctile_label , &
SOURCE_damage_anisoBrittle_label , &
SOURCE_damage_anisoDuctile_label , &
2014-05-08 20:25:19 +05:30
plasticState , &
2015-10-14 00:22:01 +05:30
sourceState
2014-08-08 16:34:40 +05:30
2014-12-08 21:25:30 +05:30
use plastic_none
2016-01-09 01:15:20 +05:30
use plastic_isotropic
2014-12-08 21:25:30 +05:30
use plastic_phenopowerlaw
2017-09-30 04:06:28 +05:30
use plastic_kinehardening
2014-12-08 21:25:30 +05:30
use plastic_dislotwin
2015-01-15 16:26:15 +05:30
use plastic_disloucla
2014-12-08 21:25:30 +05:30
use plastic_nonlocal
2015-05-28 22:32:23 +05:30
use source_thermal_dissipation
2015-07-27 16:39:37 +05:30
use source_thermal_externalheat
2015-05-28 22:32:23 +05:30
use source_damage_isoBrittle
use source_damage_isoDuctile
use source_damage_anisoBrittle
use source_damage_anisoDuctile
use kinematics_cleavage_opening
use kinematics_slipplane_opening
use kinematics_thermal_expansion
2014-10-11 02:25:09 +05:30
2013-01-16 15:44:57 +05:30
implicit none
2018-08-21 01:41:53 +05:30
integer ( pInt ) , parameter :: FILEUNIT = 204_pInt
2013-01-16 15:44:57 +05:30
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
o , & !< counter in output loop
2018-06-30 04:15:04 +05:30
ph , & !< counter in phase loop
2016-01-16 22:57:19 +05:30
s , & !< counter in source loop
ins !< instance of plasticity/source
2014-05-08 20:25:19 +05:30
2013-01-16 15:44:57 +05:30
integer ( pInt ) , dimension ( : , : ) , pointer :: thisSize
character ( len = 64 ) , dimension ( : , : ) , pointer :: thisOutput
2013-11-27 13:34:05 +05:30
character ( len = 32 ) :: outputName !< name of output, intermediate fix until HDF5 output is ready
2015-05-28 22:32:23 +05:30
logical :: knownPlasticity , knownSource , nonlocalConstitutionPresent
2013-05-24 19:13:44 +05:30
nonlocalConstitutionPresent = . false .
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2016-01-16 22:57:19 +05:30
! open material.config
2013-12-12 22:39:59 +05:30
if ( . not . IO_open_jobFile_stat ( FILEUNIT , material_localFileExt ) ) & ! no local material configuration present...
call IO_open_file ( FILEUNIT , material_configFile ) ! ... open material.config file
2016-01-16 22:57:19 +05:30
!--------------------------------------------------------------------------------------------------
! parse plasticities from config file
2014-12-08 21:25:30 +05:30
if ( any ( phase_plasticity == PLASTICITY_NONE_ID ) ) call plastic_none_init
2018-06-01 13:54:42 +05:30
if ( any ( phase_plasticity == PLASTICITY_ISOTROPIC_ID ) ) call plastic_isotropic_init
2018-06-02 02:26:20 +05:30
if ( any ( phase_plasticity == PLASTICITY_PHENOPOWERLAW_ID ) ) call plastic_phenopowerlaw_init
2018-12-13 14:21:43 +05:30
if ( any ( phase_plasticity == PLASTICITY_KINEHARDENING_ID ) ) call plastic_kinehardening_init
2018-12-11 06:17:13 +05:30
if ( any ( phase_plasticity == PLASTICITY_DISLOTWIN_ID ) ) call plastic_dislotwin_init
2018-12-05 02:03:32 +05:30
if ( any ( phase_plasticity == PLASTICITY_DISLOUCLA_ID ) ) call plastic_disloucla_init
2014-07-08 20:28:23 +05:30
if ( any ( phase_plasticity == PLASTICITY_NONLOCAL_ID ) ) then
2014-12-08 21:25:30 +05:30
call plastic_nonlocal_init ( FILEUNIT )
call plastic_nonlocal_stateInit ( )
2014-07-08 20:28:23 +05:30
endif
2014-09-23 16:08:20 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
! parse source mechanisms from config file
2016-01-16 22:57:19 +05:30
call IO_checkAndRewind ( FILEUNIT )
2015-05-28 22:32:23 +05:30
if ( any ( phase_source == SOURCE_thermal_dissipation_ID ) ) call source_thermal_dissipation_init ( FILEUNIT )
2015-07-27 16:39:37 +05:30
if ( any ( phase_source == SOURCE_thermal_externalheat_ID ) ) call source_thermal_externalheat_init ( FILEUNIT )
2015-05-28 22:32:23 +05:30
if ( any ( phase_source == SOURCE_damage_isoBrittle_ID ) ) call source_damage_isoBrittle_init ( FILEUNIT )
if ( any ( phase_source == SOURCE_damage_isoDuctile_ID ) ) call source_damage_isoDuctile_init ( FILEUNIT )
if ( any ( phase_source == SOURCE_damage_anisoBrittle_ID ) ) call source_damage_anisoBrittle_init ( FILEUNIT )
if ( any ( phase_source == SOURCE_damage_anisoDuctile_ID ) ) call source_damage_anisoDuctile_init ( FILEUNIT )
2015-10-14 00:22:01 +05:30
2014-09-23 16:08:20 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
! parse kinematic mechanisms from config file
2016-01-16 22:57:19 +05:30
call IO_checkAndRewind ( FILEUNIT )
2015-05-28 22:32:23 +05:30
if ( any ( phase_kinematics == KINEMATICS_cleavage_opening_ID ) ) call kinematics_cleavage_opening_init ( FILEUNIT )
if ( any ( phase_kinematics == KINEMATICS_slipplane_opening_ID ) ) call kinematics_slipplane_opening_init ( FILEUNIT )
if ( any ( phase_kinematics == KINEMATICS_thermal_expansion_ID ) ) call kinematics_thermal_expansion_init ( FILEUNIT )
2014-09-23 16:08:20 +05:30
close ( FILEUNIT )
2018-06-27 00:01:52 +05:30
call config_deallocate ( 'material.config/phase' )
2017-11-21 19:38:45 +05:30
write ( 6 , '(/,a)' ) ' <<<+- constitutive init -+>>>'
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-10-09 18:04:57 +05:30
#include "compilation_info.f90"
2015-10-14 00:22:01 +05:30
2017-11-21 19:38:45 +05:30
mainProcess : if ( worldrank == 0 ) then
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2016-01-16 22:57:19 +05:30
! write description file for constitutive output
2015-10-14 00:22:01 +05:30
call IO_write_jobFile ( FILEUNIT , 'outputConstitutive' )
2018-06-30 04:15:04 +05:30
PhaseLoop : do ph = 1_pInt , material_Nphase
activePhase : if ( any ( material_phase == ph ) ) then
ins = phase_plasticityInstance ( ph )
2016-01-16 22:57:19 +05:30
knownPlasticity = . true . ! assume valid
2018-06-30 04:15:04 +05:30
plasticityType : select case ( phase_plasticity ( ph ) )
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONE_ID ) plasticityType
2015-04-21 21:41:30 +05:30
outputName = PLASTICITY_NONE_label
2016-01-16 22:57:19 +05:30
thisOutput = > null ( )
2016-04-13 18:23:13 +05:30
thisSize = > null ( )
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_ISOTROPIC_ID ) plasticityType
2016-01-09 01:15:20 +05:30
outputName = PLASTICITY_ISOTROPIC_label
thisOutput = > plastic_isotropic_output
thisSize = > plastic_isotropic_sizePostResult
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_PHENOPOWERLAW_ID ) plasticityType
2015-04-21 21:41:30 +05:30
outputName = PLASTICITY_PHENOPOWERLAW_label
thisOutput = > plastic_phenopowerlaw_output
thisSize = > plastic_phenopowerlaw_sizePostResult
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_KINEHARDENING_ID ) plasticityType
outputName = PLASTICITY_KINEHARDENING_label
thisOutput = > plastic_kinehardening_output
thisSize = > plastic_kinehardening_sizePostResult
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2015-04-21 21:41:30 +05:30
outputName = PLASTICITY_DISLOTWIN_label
thisOutput = > plastic_dislotwin_output
thisSize = > plastic_dislotwin_sizePostResult
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOUCLA_ID ) plasticityType
2015-04-21 21:41:30 +05:30
outputName = PLASTICITY_DISLOUCLA_label
thisOutput = > plastic_disloucla_output
2015-10-14 00:22:01 +05:30
thisSize = > plastic_disloucla_sizePostResult
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2015-04-21 21:41:30 +05:30
outputName = PLASTICITY_NONLOCAL_label
thisOutput = > plastic_nonlocal_output
thisSize = > plastic_nonlocal_sizePostResult
2016-01-16 22:57:19 +05:30
case default plasticityType
2015-04-21 21:41:30 +05:30
knownPlasticity = . false .
2016-01-16 22:57:19 +05:30
end select plasticityType
2018-06-30 04:15:04 +05:30
write ( FILEUNIT , '(/,a,/)' ) '[' / / trim ( phase_name ( ph ) ) / / ']'
2015-04-21 21:41:30 +05:30
if ( knownPlasticity ) then
write ( FILEUNIT , '(a)' ) '(plasticity)' / / char ( 9 ) / / trim ( outputName )
2018-06-30 04:15:04 +05:30
if ( phase_plasticity ( ph ) / = PLASTICITY_NONE_ID ) then
2018-06-02 00:24:40 +05:30
OutputPlasticityLoop : do o = 1_pInt , size ( thisOutput ( : , ins ) )
2018-06-02 02:44:43 +05:30
if ( len ( trim ( thisOutput ( o , ins ) ) ) > 0_pInt ) &
write ( FILEUNIT , '(a,i4)' ) trim ( thisOutput ( o , ins ) ) / / char ( 9 ) , thisSize ( o , ins )
2016-01-16 22:57:19 +05:30
enddo OutputPlasticityLoop
2015-04-21 21:41:30 +05:30
endif
2015-03-25 21:32:30 +05:30
endif
2018-06-30 04:15:04 +05:30
SourceLoop : do s = 1_pInt , phase_Nsources ( ph )
2016-01-16 22:57:19 +05:30
knownSource = . true . ! assume valid
2018-06-30 04:15:04 +05:30
sourceType : select case ( phase_source ( s , ph ) )
2016-01-16 22:57:19 +05:30
case ( SOURCE_thermal_dissipation_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_thermal_dissipation_instance ( ph )
2015-05-28 22:32:23 +05:30
outputName = SOURCE_thermal_dissipation_label
thisOutput = > source_thermal_dissipation_output
thisSize = > source_thermal_dissipation_sizePostResult
2016-01-16 22:57:19 +05:30
case ( SOURCE_thermal_externalheat_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_thermal_externalheat_instance ( ph )
2015-07-27 16:39:37 +05:30
outputName = SOURCE_thermal_externalheat_label
thisOutput = > source_thermal_externalheat_output
thisSize = > source_thermal_externalheat_sizePostResult
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_isoBrittle_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_damage_isoBrittle_instance ( ph )
2015-05-28 22:32:23 +05:30
outputName = SOURCE_damage_isoBrittle_label
thisOutput = > source_damage_isoBrittle_output
thisSize = > source_damage_isoBrittle_sizePostResult
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_isoDuctile_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_damage_isoDuctile_instance ( ph )
2015-05-28 22:32:23 +05:30
outputName = SOURCE_damage_isoDuctile_label
thisOutput = > source_damage_isoDuctile_output
thisSize = > source_damage_isoDuctile_sizePostResult
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoBrittle_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_damage_anisoBrittle_instance ( ph )
2015-05-28 22:32:23 +05:30
outputName = SOURCE_damage_anisoBrittle_label
thisOutput = > source_damage_anisoBrittle_output
thisSize = > source_damage_anisoBrittle_sizePostResult
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoDuctile_ID ) sourceType
2018-06-30 04:15:04 +05:30
ins = source_damage_anisoDuctile_instance ( ph )
2015-05-28 22:32:23 +05:30
outputName = SOURCE_damage_anisoDuctile_label
thisOutput = > source_damage_anisoDuctile_output
thisSize = > source_damage_anisoDuctile_sizePostResult
2016-01-16 22:57:19 +05:30
case default sourceType
2015-05-28 22:32:23 +05:30
knownSource = . false .
2016-01-16 22:57:19 +05:30
end select sourceType
2015-05-28 22:32:23 +05:30
if ( knownSource ) then
write ( FILEUNIT , '(a)' ) '(source)' / / char ( 9 ) / / trim ( outputName )
2018-06-02 00:24:40 +05:30
OutputSourceLoop : do o = 1_pInt , size ( thisOutput ( : , ins ) )
2018-06-02 02:44:43 +05:30
if ( len ( trim ( thisOutput ( o , ins ) ) ) > 0_pInt ) &
write ( FILEUNIT , '(a,i4)' ) trim ( thisOutput ( o , ins ) ) / / char ( 9 ) , thisSize ( o , ins )
2016-01-16 22:57:19 +05:30
enddo OutputSourceLoop
2015-04-21 21:41:30 +05:30
endif
2016-01-16 22:57:19 +05:30
enddo SourceLoop
endif activePhase
enddo PhaseLoop
2015-03-25 21:32:30 +05:30
close ( FILEUNIT )
2016-01-16 22:57:19 +05:30
endif mainProcess
2015-10-14 00:22:01 +05:30
2015-05-28 22:32:23 +05:30
constitutive_plasticity_maxSizeDotState = 0_pInt
constitutive_plasticity_maxSizePostResults = 0_pInt
constitutive_source_maxSizeDotState = 0_pInt
constitutive_source_maxSizePostResults = 0_pInt
2014-07-03 18:47:29 +05:30
2018-06-30 04:15:04 +05:30
PhaseLoop2 : do ph = 1_pInt , material_Nphase
2016-01-16 22:57:19 +05:30
!--------------------------------------------------------------------------------------------------
! partition and inititalize state
2018-06-30 04:15:04 +05:30
plasticState ( ph ) % partionedState0 = plasticState ( ph ) % state0
plasticState ( ph ) % state = plasticState ( ph ) % partionedState0
forall ( s = 1_pInt : phase_Nsources ( ph ) )
sourceState ( ph ) % p ( s ) % partionedState0 = sourceState ( ph ) % p ( s ) % state0
sourceState ( ph ) % p ( s ) % state = sourceState ( ph ) % p ( s ) % partionedState0
2016-01-16 22:57:19 +05:30
end forall
!--------------------------------------------------------------------------------------------------
! determine max size of state and output
2015-05-28 22:32:23 +05:30
constitutive_plasticity_maxSizeDotState = max ( constitutive_plasticity_maxSizeDotState , &
2018-06-30 04:15:04 +05:30
plasticState ( ph ) % sizeDotState )
2015-05-28 22:32:23 +05:30
constitutive_plasticity_maxSizePostResults = max ( constitutive_plasticity_maxSizePostResults , &
2018-06-30 04:15:04 +05:30
plasticState ( ph ) % sizePostResults )
2015-05-28 22:32:23 +05:30
constitutive_source_maxSizeDotState = max ( constitutive_source_maxSizeDotState , &
2018-06-30 04:15:04 +05:30
maxval ( sourceState ( ph ) % p ( : ) % sizeDotState ) )
2015-05-28 22:32:23 +05:30
constitutive_source_maxSizePostResults = max ( constitutive_source_maxSizePostResults , &
2018-06-30 04:15:04 +05:30
maxval ( sourceState ( ph ) % p ( : ) % sizePostResults ) )
2014-07-03 18:47:29 +05:30
enddo PhaseLoop2
2014-06-24 14:54:59 +05:30
2014-05-08 20:25:19 +05:30
2012-10-02 18:23:25 +05:30
end subroutine constitutive_init
2009-03-06 15:32:36 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns the homogenize elasticity matrix
2018-08-25 17:25:39 +05:30
!> ToDo: homogenizedC66 would be more consistent
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2014-07-02 17:57:39 +05:30
function constitutive_homogenizedC ( ipc , ip , el )
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2013-01-16 15:44:57 +05:30
use material , only : &
phase_plasticity , &
2013-11-27 13:34:05 +05:30
material_phase , &
2014-08-08 16:34:40 +05:30
PLASTICITY_DISLOTWIN_ID , &
2015-04-21 21:41:30 +05:30
PLASTICITY_DISLOUCLA_ID
2014-12-08 21:25:30 +05:30
use plastic_dislotwin , only : &
plastic_dislotwin_homogenizedC
2014-03-09 02:20:31 +05:30
use lattice , only : &
lattice_C66
2009-03-04 19:31:36 +05:30
implicit none
real ( pReal ) , dimension ( 6 , 6 ) :: constitutive_homogenizedC
2013-01-16 15:44:57 +05:30
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2009-03-06 15:32:36 +05:30
2016-01-16 22:57:19 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2014-12-08 21:25:30 +05:30
constitutive_homogenizedC = plastic_dislotwin_homogenizedC ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
case default plasticityType
2014-07-02 17:57:39 +05:30
constitutive_homogenizedC = lattice_C66 ( 1 : 6 , 1 : 6 , material_phase ( ipc , ip , el ) )
2016-01-16 22:57:19 +05:30
end select plasticityType
2009-03-06 15:32:36 +05:30
2012-10-02 18:23:25 +05:30
end function constitutive_homogenizedC
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calls microstructure function of the different constitutive models
!--------------------------------------------------------------------------------------------------
2017-02-04 01:05:05 +05:30
subroutine constitutive_microstructure ( orientations , Fe , Fp , ipc , ip , el )
2014-07-02 17:57:39 +05:30
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2013-01-16 15:44:57 +05:30
use material , only : &
2018-12-05 02:23:22 +05:30
phasememberAt , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
2018-12-05 02:23:22 +05:30
phase_plasticityInstance , &
2013-11-27 13:34:05 +05:30
material_phase , &
2018-10-04 10:09:03 +05:30
material_homogenizationAt , &
2015-05-28 22:32:23 +05:30
temperature , &
thermalMapping , &
2014-10-10 18:12:12 +05:30
PLASTICITY_dislotwin_ID , &
2015-01-15 16:26:15 +05:30
PLASTICITY_disloucla_ID , &
2017-10-08 17:11:54 +05:30
PLASTICITY_nonlocal_ID
2014-12-08 21:25:30 +05:30
use plastic_nonlocal , only : &
plastic_nonlocal_microstructure
use plastic_dislotwin , only : &
2019-01-27 13:05:07 +05:30
plastic_dislotwin_dependentState
2018-12-05 02:23:22 +05:30
use plastic_disloUCLA , only : &
plastic_disloUCLA_dependentState
2014-06-03 19:16:42 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
Fe , & !< elastic deformation gradient
2015-05-28 22:32:23 +05:30
Fp !< plastic deformation gradient
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
ho , & !< homogenization
2018-12-05 02:23:22 +05:30
tme , & !< thermal member position
instance , of
2015-10-14 00:22:01 +05:30
real ( pReal ) , intent ( in ) , dimension ( : , : , : , : ) :: &
2017-02-04 01:05:05 +05:30
orientations !< crystal orientations as quaternions
2016-01-16 22:57:19 +05:30
2018-10-04 10:09:03 +05:30
ho = material_homogenizationAt ( el )
2016-01-16 22:57:19 +05:30
tme = thermalMapping ( ho ) % p ( ip , el )
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2019-01-27 13:05:07 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_dislotwin_dependentState ( temperature ( ho ) % p ( tme ) , instance , of )
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOUCLA_ID ) plasticityType
2018-12-05 02:23:22 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2018-12-05 03:00:07 +05:30
call plastic_disloUCLA_dependentState ( instance , of )
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2015-05-28 22:32:23 +05:30
call plastic_nonlocal_microstructure ( Fe , Fp , ip , el )
2016-01-16 22:57:19 +05:30
end select plasticityType
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
end subroutine constitutive_microstructure
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief contains the constitutive equation for calculating the velocity gradient
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-28 18:48:55 +05:30
subroutine constitutive_LpAndItsTangents ( Lp , dLp_dS , dLp_dFi , S6 , Fi , ipc , ip , el )
2014-07-02 17:57:39 +05:30
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2013-10-16 18:34:59 +05:30
use math , only : &
2015-03-06 18:39:00 +05:30
math_mul33x33 , &
2019-01-27 15:18:27 +05:30
math_6toSym33 , &
math_sym33to6 , &
math_99to3333
2013-01-16 15:44:57 +05:30
use material , only : &
2018-09-17 01:46:06 +05:30
phasememberAt , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
2018-09-17 01:46:06 +05:30
phase_plasticityInstance , &
2013-11-27 13:34:05 +05:30
material_phase , &
2018-10-04 10:09:03 +05:30
material_homogenizationAt , &
2015-05-28 22:32:23 +05:30
temperature , &
thermalMapping , &
2013-11-27 13:34:05 +05:30
PLASTICITY_NONE_ID , &
2016-01-09 01:15:20 +05:30
PLASTICITY_ISOTROPIC_ID , &
2013-11-27 13:34:05 +05:30
PLASTICITY_PHENOPOWERLAW_ID , &
2017-09-30 04:06:28 +05:30
PLASTICITY_KINEHARDENING_ID , &
2013-11-27 13:34:05 +05:30
PLASTICITY_DISLOTWIN_ID , &
2015-01-15 16:26:15 +05:30
PLASTICITY_DISLOUCLA_ID , &
2013-11-27 13:34:05 +05:30
PLASTICITY_NONLOCAL_ID
2016-01-09 01:15:20 +05:30
use plastic_isotropic , only : &
plastic_isotropic_LpAndItsTangent
2014-12-08 21:25:30 +05:30
use plastic_phenopowerlaw , only : &
2015-01-29 19:26:09 +05:30
plastic_phenopowerlaw_LpAndItsTangent
2017-09-30 04:06:28 +05:30
use plastic_kinehardening , only : &
plastic_kinehardening_LpAndItsTangent
2014-12-08 21:25:30 +05:30
use plastic_dislotwin , only : &
2015-01-29 19:26:09 +05:30
plastic_dislotwin_LpAndItsTangent
2015-01-15 16:26:15 +05:30
use plastic_disloucla , only : &
2015-01-29 19:26:09 +05:30
plastic_disloucla_LpAndItsTangent
2014-12-08 21:25:30 +05:30
use plastic_nonlocal , only : &
2015-01-29 19:26:09 +05:30
plastic_nonlocal_LpAndItsTangent
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 6 ) :: &
2018-08-28 18:32:58 +05:30
S6 !< 2nd Piola-Kirchhoff stress (vector notation)
2015-03-06 18:39:00 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
Fi !< intermediate deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 ) :: &
Lp !< plastic velocity gradient
2015-03-06 18:39:00 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , 3 , 3 ) :: &
2018-08-28 18:48:55 +05:30
dLp_dS , &
2018-08-25 17:51:06 +05:30
dLp_dFi !< derivative of Lp with respect to Fi
2018-08-28 18:48:55 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
dLp_dMp !< derivative of Lp with respect to Mandel stress
2015-03-06 18:39:00 +05:30
real ( pReal ) , dimension ( 9 , 9 ) :: &
2018-08-28 18:48:55 +05:30
dLp_dMp99 !< derivative of Lp with respect to Mstar (matrix notation)
2015-03-06 18:39:00 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
2018-08-28 18:48:55 +05:30
Mp , & !< Mandel stress work conjugate with Lp
2018-08-28 18:32:58 +05:30
S !< 2nd Piola-Kirchhoff stress
2015-03-06 18:39:00 +05:30
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
ho , & !< homogenization
tme !< thermal member position
integer ( pInt ) :: &
2018-09-17 01:46:06 +05:30
i , j , instance , of
2016-01-16 22:57:19 +05:30
2018-10-04 10:09:03 +05:30
ho = material_homogenizationAt ( el )
2016-01-16 22:57:19 +05:30
tme = thermalMapping ( ho ) % p ( ip , el )
2014-10-28 08:12:25 +05:30
2019-01-27 15:18:27 +05:30
S = math_6toSym33 ( S6 )
2018-08-28 18:48:55 +05:30
Mp = math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , Fi ) , S )
2015-10-14 00:22:01 +05:30
2016-01-16 22:57:19 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONE_ID ) plasticityType
2013-10-16 18:34:59 +05:30
Lp = 0.0_pReal
2018-08-28 18:48:55 +05:30
dLp_dMp = 0.0_pReal
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_ISOTROPIC_ID ) plasticityType
2018-12-30 20:39:51 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_isotropic_LpAndItsTangent ( Lp , dLp_dMp , Mp , instance , of )
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_PHENOPOWERLAW_ID ) plasticityType
2018-09-17 01:46:06 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2018-10-02 01:44:54 +05:30
call plastic_phenopowerlaw_LpAndItsTangent ( Lp , dLp_dMp , Mp , instance , of )
2018-08-25 17:51:06 +05:30
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_KINEHARDENING_ID ) plasticityType
2018-12-13 14:04:40 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_kinehardening_LpAndItsTangent ( Lp , dLp_dMp , Mp , instance , of )
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2019-01-27 15:18:27 +05:30
call plastic_nonlocal_LpAndItsTangent ( Lp , dLp_dMp99 , math_sym33to6 ( Mp ) , &
2017-09-30 04:06:28 +05:30
temperature ( ho ) % p ( tme ) , ip , el )
2019-01-27 15:18:27 +05:30
dLp_dMp = math_99to3333 ( dLp_dMp99 ) ! ToDo: We revert here the last statement in plastic_xx_LpAndItsTanget
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2018-10-02 01:44:54 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_dislotwin_LpAndItsTangent ( Lp , dLp_dMp , Mp , temperature ( ho ) % p ( tme ) , instance , of )
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOUCLA_ID ) plasticityType
2018-12-05 02:03:32 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_disloucla_LpAndItsTangent ( Lp , dLp_dMp , Mp , temperature ( ho ) % p ( tme ) , instance , of )
2018-08-25 17:51:06 +05:30
2016-01-16 22:57:19 +05:30
end select plasticityType
2015-03-06 18:39:00 +05:30
2018-10-09 12:05:54 +05:30
#ifdef __INTEL_COMPILER
forall ( i = 1_pInt : 3_pInt , j = 1_pInt : 3_pInt )
#else
2018-08-25 18:06:21 +05:30
do concurrent ( i = 1_pInt : 3_pInt , j = 1_pInt : 3_pInt )
2018-10-09 12:05:54 +05:30
#endif
2018-08-28 18:48:55 +05:30
dLp_dFi ( i , j , 1 : 3 , 1 : 3 ) = math_mul33x33 ( math_mul33x33 ( Fi , S ) , transpose ( dLp_dMp ( i , j , 1 : 3 , 1 : 3 ) ) ) + &
math_mul33x33 ( math_mul33x33 ( Fi , dLp_dMp ( i , j , 1 : 3 , 1 : 3 ) ) , S )
dLp_dS ( i , j , 1 : 3 , 1 : 3 ) = math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , Fi ) , dLp_dMp ( i , j , 1 : 3 , 1 : 3 ) ) ! ToDo: @PS: why not: dLp_dMp:(FiT Fi)
2018-10-09 12:05:54 +05:30
#ifdef __INTEL_COMPILER
end forall
#else
enddo
#endif
2015-10-14 00:22:01 +05:30
2018-08-28 18:37:39 +05:30
end subroutine constitutive_LpAndItsTangents
2009-03-06 15:32:36 +05:30
2014-11-01 00:33:08 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief contains the constitutive equation for calculating the velocity gradient
2018-12-30 17:05:26 +05:30
! ToDo: MD: S is Mi?
2014-11-01 00:33:08 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-28 18:37:39 +05:30
subroutine constitutive_LiAndItsTangents ( Li , dLi_dS , dLi_dFi , S6 , Fi , ipc , ip , el )
2014-11-01 00:33:08 +05:30
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2015-03-06 18:39:00 +05:30
use math , only : &
math_I3 , &
math_inv33 , &
math_det33 , &
2018-12-30 17:05:26 +05:30
math_mul33x33 , &
2019-01-27 15:18:27 +05:30
math_6toSym33
2014-11-01 00:33:08 +05:30
use material , only : &
2018-12-30 18:31:05 +05:30
phasememberAt , &
phase_plasticity , &
phase_plasticityInstance , &
2016-01-09 01:15:20 +05:30
phase_plasticity , &
material_phase , &
2015-05-28 22:32:23 +05:30
phase_kinematics , &
phase_Nkinematics , &
2016-01-09 01:15:20 +05:30
PLASTICITY_isotropic_ID , &
2015-05-28 22:32:23 +05:30
KINEMATICS_cleavage_opening_ID , &
KINEMATICS_slipplane_opening_ID , &
2018-12-22 13:07:58 +05:30
KINEMATICS_thermal_expansion_ID
2016-01-09 01:15:20 +05:30
use plastic_isotropic , only : &
plastic_isotropic_LiAndItsTangent
2015-05-28 22:32:23 +05:30
use kinematics_cleavage_opening , only : &
kinematics_cleavage_opening_LiAndItsTangent
use kinematics_slipplane_opening , only : &
kinematics_slipplane_opening_LiAndItsTangent
use kinematics_thermal_expansion , only : &
kinematics_thermal_expansion_LiAndItsTangent
2015-10-14 00:22:01 +05:30
2014-11-01 00:33:08 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2014-11-01 00:33:08 +05:30
real ( pReal ) , intent ( in ) , dimension ( 6 ) :: &
2018-08-28 18:32:58 +05:30
S6 !< 2nd Piola-Kirchhoff stress (vector notation)
2014-11-01 00:33:08 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
2015-05-28 22:32:23 +05:30
Fi !< intermediate deformation gradient
2014-11-01 00:33:08 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 ) :: &
Li !< intermediate velocity gradient
2015-03-06 18:39:00 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , 3 , 3 ) :: &
2018-08-28 18:32:58 +05:30
dLi_dS , & !< derivative of Li with respect to S
2018-08-25 18:12:44 +05:30
dLi_dFi
2018-12-30 18:31:05 +05:30
2014-11-01 00:33:08 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
2018-12-30 18:31:05 +05:30
my_Li , & !< intermediate velocity gradient
2015-03-06 18:39:00 +05:30
FiInv , &
temp_33
2018-12-30 18:31:05 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
my_dLi_dS
2015-03-06 18:39:00 +05:30
real ( pReal ) :: &
detFi
integer ( pInt ) :: &
2018-12-30 18:31:05 +05:30
k , i , j , &
instance , of
2014-11-01 00:33:08 +05:30
Li = 0.0_pReal
2018-08-28 18:32:58 +05:30
dLi_dS = 0.0_pReal
2018-08-28 18:48:55 +05:30
dLi_dFi = 0.0_pReal
2015-10-14 00:22:01 +05:30
2016-01-16 22:57:19 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
case ( PLASTICITY_isotropic_ID ) plasticityType
2018-12-30 18:31:05 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2019-01-27 15:18:27 +05:30
call plastic_isotropic_LiAndItsTangent ( my_Li , my_dLi_dS , math_6toSym33 ( S6 ) , instance , of )
2016-01-16 22:57:19 +05:30
case default plasticityType
2016-01-09 01:15:20 +05:30
my_Li = 0.0_pReal
2018-08-28 18:32:58 +05:30
my_dLi_dS = 0.0_pReal
2016-01-16 22:57:19 +05:30
end select plasticityType
2016-01-09 01:15:20 +05:30
Li = Li + my_Li
2018-08-28 18:32:58 +05:30
dLi_dS = dLi_dS + my_dLi_dS
2016-01-09 01:15:20 +05:30
2016-01-16 22:57:19 +05:30
KinematicsLoop : do k = 1_pInt , phase_Nkinematics ( material_phase ( ipc , ip , el ) )
kinematicsType : select case ( phase_kinematics ( k , material_phase ( ipc , ip , el ) ) )
case ( KINEMATICS_cleavage_opening_ID ) kinematicsType
2018-08-28 18:32:58 +05:30
call kinematics_cleavage_opening_LiAndItsTangent ( my_Li , my_dLi_dS , S6 , ipc , ip , el )
2016-01-16 22:57:19 +05:30
case ( KINEMATICS_slipplane_opening_ID ) kinematicsType
2018-08-28 18:32:58 +05:30
call kinematics_slipplane_opening_LiAndItsTangent ( my_Li , my_dLi_dS , S6 , ipc , ip , el )
2016-01-16 22:57:19 +05:30
case ( KINEMATICS_thermal_expansion_ID ) kinematicsType
2018-08-28 18:32:58 +05:30
call kinematics_thermal_expansion_LiAndItsTangent ( my_Li , my_dLi_dS , ipc , ip , el )
2016-01-16 22:57:19 +05:30
case default kinematicsType
2015-05-28 22:32:23 +05:30
my_Li = 0.0_pReal
2018-08-28 18:32:58 +05:30
my_dLi_dS = 0.0_pReal
2016-01-16 22:57:19 +05:30
end select kinematicsType
2015-05-28 22:32:23 +05:30
Li = Li + my_Li
2018-08-28 18:32:58 +05:30
dLi_dS = dLi_dS + my_dLi_dS
2016-01-16 22:57:19 +05:30
enddo KinematicsLoop
2014-11-01 00:33:08 +05:30
2015-03-06 18:39:00 +05:30
FiInv = math_inv33 ( Fi )
detFi = math_det33 ( Fi )
Li = math_mul33x33 ( math_mul33x33 ( Fi , Li ) , FiInv ) * detFi !< push forward to intermediate configuration
temp_33 = math_mul33x33 ( FiInv , Li )
2018-10-13 12:05:36 +05:30
do i = 1_pInt , 3_pInt ; do j = 1_pInt , 3_pInt
2018-08-28 18:48:55 +05:30
dLi_dS ( 1 : 3 , 1 : 3 , i , j ) = math_mul33x33 ( math_mul33x33 ( Fi , dLi_dS ( 1 : 3 , 1 : 3 , i , j ) ) , FiInv ) * detFi
dLi_dFi ( 1 : 3 , 1 : 3 , i , j ) = dLi_dFi ( 1 : 3 , 1 : 3 , i , j ) + Li * FiInv ( j , i )
dLi_dFi ( 1 : 3 , i , 1 : 3 , j ) = dLi_dFi ( 1 : 3 , i , 1 : 3 , j ) + math_I3 * temp_33 ( j , i ) + Li * FiInv ( j , i )
2018-10-13 12:05:36 +05:30
end do ; end do
2015-10-14 00:22:01 +05:30
2018-08-28 18:37:39 +05:30
end subroutine constitutive_LiAndItsTangents
2014-11-01 00:33:08 +05:30
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
2015-07-24 20:17:18 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief collects initial intermediate deformation gradient
2015-07-24 20:17:18 +05:30
!--------------------------------------------------------------------------------------------------
pure function constitutive_initialFi ( ipc , ip , el )
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2015-07-24 20:17:18 +05:30
use math , only : &
math_I3 , &
math_inv33 , &
math_mul33x33
use material , only : &
phase_kinematics , &
phase_Nkinematics , &
material_phase , &
2018-12-22 13:07:58 +05:30
KINEMATICS_thermal_expansion_ID
2015-07-24 20:17:18 +05:30
use kinematics_thermal_expansion , only : &
kinematics_thermal_expansion_initialStrain
2015-10-14 00:22:01 +05:30
2015-07-24 20:17:18 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2015-07-24 20:17:18 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
constitutive_initialFi !< composite initial intermediate deformation gradient
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
k !< counter in kinematics loop
2015-07-24 20:17:18 +05:30
constitutive_initialFi = math_I3
2015-10-14 00:22:01 +05:30
2016-01-16 22:57:19 +05:30
KinematicsLoop : do k = 1_pInt , phase_Nkinematics ( material_phase ( ipc , ip , el ) ) !< Warning: small initial strain assumption
kinematicsType : select case ( phase_kinematics ( k , material_phase ( ipc , ip , el ) ) )
case ( KINEMATICS_thermal_expansion_ID ) kinematicsType
2015-07-24 20:17:18 +05:30
constitutive_initialFi = &
constitutive_initialFi + kinematics_thermal_expansion_initialStrain ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
end select kinematicsType
enddo KinematicsLoop
2015-10-14 00:22:01 +05:30
2015-07-24 20:17:18 +05:30
end function constitutive_initialFi
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief returns the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to
2018-08-28 18:24:36 +05:30
!> the elastic/intermediate deformation gradients depending on the selected elastic law
!! (so far no case switch because only Hooke is implemented)
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-28 18:24:36 +05:30
subroutine constitutive_SandItsTangents ( S , dS_dFe , dS_dFi , Fe , Fi , ipc , ip , el )
2014-07-02 17:57:39 +05:30
use prec , only : &
pReal
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
2015-03-06 18:39:00 +05:30
Fe , & !< elastic deformation gradient
Fi !< intermediate deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 ) :: &
2018-08-28 18:24:36 +05:30
S !< 2nd Piola-Kirchhoff stress tensor
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , 3 , 3 ) :: &
2018-08-28 18:24:36 +05:30
dS_dFe , & !< derivative of 2nd P-K stress with respect to elastic deformation gradient
dS_dFi !< derivative of 2nd P-K stress with respect to intermediate deformation gradient
2015-10-14 00:22:01 +05:30
2018-08-28 18:24:36 +05:30
call constitutive_hooke_SandItsTangents ( S , dS_dFe , dS_dFi , Fe , Fi , ipc , ip , el )
2013-10-14 16:24:45 +05:30
2015-10-14 00:22:01 +05:30
2018-08-28 18:24:36 +05:30
end subroutine constitutive_SandItsTangents
2012-03-15 15:21:33 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief returns the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to
2018-08-28 18:24:36 +05:30
!> the elastic and intermeidate deformation gradients using Hookes law
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-28 18:24:36 +05:30
subroutine constitutive_hooke_SandItsTangents ( S , dS_dFe , dS_dFi , Fe , Fi , ipc , ip , el )
2014-07-02 17:57:39 +05:30
use prec , only : &
pReal
use math , only : &
math_mul33x33 , &
math_mul3333xx33 , &
2019-01-27 15:18:27 +05:30
math_66toSym3333 , &
2014-08-09 03:59:38 +05:30
math_I3
2015-05-28 22:32:23 +05:30
use material , only : &
material_phase , &
2018-10-04 10:09:03 +05:30
material_homogenizationAt , &
2015-05-28 22:32:23 +05:30
phase_NstiffnessDegradations , &
phase_stiffnessDegradation , &
damage , &
damageMapping , &
2018-12-30 15:11:11 +05:30
STIFFNESS_DEGRADATION_damage_ID
2012-11-07 21:13:29 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
2015-03-06 18:39:00 +05:30
Fe , & !< elastic deformation gradient
Fi !< intermediate deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 ) :: &
2018-08-28 18:24:36 +05:30
S !< 2nd Piola-Kirchhoff stress tensor in lattice configuration
2015-03-06 18:39:00 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , 3 , 3 ) :: &
2018-08-28 18:24:36 +05:30
dS_dFe , & !< derivative of 2nd P-K stress with respect to elastic deformation gradient
dS_dFi !< derivative of 2nd P-K stress with respect to intermediate deformation gradient
2015-03-06 18:39:00 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: E
2014-08-09 03:59:38 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: C
2016-01-16 22:57:19 +05:30
integer ( pInt ) :: &
ho , & !< homogenization
d !< counter in degradation loop
integer ( pInt ) :: &
i , j
2018-10-04 10:09:03 +05:30
ho = material_homogenizationAt ( el )
2019-01-27 15:18:27 +05:30
C = math_66toSym3333 ( constitutive_homogenizedC ( ipc , ip , el ) )
2016-01-16 22:57:19 +05:30
DegradationLoop : do d = 1_pInt , phase_NstiffnessDegradations ( material_phase ( ipc , ip , el ) )
degradationType : select case ( phase_stiffnessDegradation ( d , material_phase ( ipc , ip , el ) ) )
case ( STIFFNESS_DEGRADATION_damage_ID ) degradationType
C = C * damage ( ho ) % p ( damageMapping ( ho ) % p ( ip , el ) ) ** 2_pInt
end select degradationType
enddo DegradationLoop
E = 0.5_pReal * ( math_mul33x33 ( transpose ( Fe ) , Fe ) - math_I3 ) !< Green-Lagrange strain in unloaded configuration
2018-08-28 18:24:36 +05:30
S = math_mul3333xx33 ( C , math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , E ) , Fi ) ) !< 2PK stress in lattice configuration in work conjugate with GL strain pulled back to lattice configuration
2015-10-14 00:22:01 +05:30
2018-08-28 18:24:36 +05:30
dS_dFe = 0.0_pReal
2015-03-06 18:39:00 +05:30
forall ( i = 1_pInt : 3_pInt , j = 1_pInt : 3_pInt )
2018-08-28 18:24:36 +05:30
dS_dFe ( i , j , 1 : 3 , 1 : 3 ) = &
math_mul33x33 ( Fe , math_mul33x33 ( math_mul33x33 ( Fi , C ( i , j , 1 : 3 , 1 : 3 ) ) , transpose ( Fi ) ) ) !< dS_ij/dFe_kl = C_ijmn * Fi_lm * Fi_on * Fe_ko
dS_dFi ( i , j , 1 : 3 , 1 : 3 ) = 2.0_pReal * math_mul33x33 ( math_mul33x33 ( E , Fi ) , C ( i , j , 1 : 3 , 1 : 3 ) ) !< dS_ij/dFi_kl = C_ijln * E_km * Fe_mn
2015-03-06 18:39:00 +05:30
end forall
2015-10-14 00:22:01 +05:30
2018-08-28 18:24:36 +05:30
end subroutine constitutive_hooke_SandItsTangents
2012-03-15 15:21:33 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-14 00:22:01 +05:30
!> @brief contains the constitutive equation for calculating the rate of change of microstructure
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2018-08-28 18:32:58 +05:30
subroutine constitutive_collectDotState ( S6 , FeArray , Fi , FpArray , subdt , subfracArray , ipc , ip , el )
2013-01-16 15:44:57 +05:30
use prec , only : &
2014-07-02 17:57:39 +05:30
pReal , &
2013-01-16 15:44:57 +05:30
pLongInt
use debug , only : &
debug_level , &
debug_constitutive , &
debug_levelBasic
2018-08-25 19:29:34 +05:30
use math , only : &
2018-09-17 02:27:50 +05:30
math_mul33x33 , &
2019-01-27 15:18:27 +05:30
math_6toSym33 , &
math_sym33to6 , &
2018-08-25 19:29:34 +05:30
math_mul33x33
2013-01-16 15:44:57 +05:30
use mesh , only : &
mesh_NcpElems , &
2013-05-08 14:53:47 +05:30
mesh_maxNips
2013-01-16 15:44:57 +05:30
use material , only : &
2018-09-17 01:46:06 +05:30
phasememberAt , &
phase_plasticityInstance , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
2015-05-28 22:32:23 +05:30
phase_source , &
phase_Nsources , &
2013-01-16 15:44:57 +05:30
material_phase , &
2018-10-04 10:09:03 +05:30
material_homogenizationAt , &
2015-05-28 22:32:23 +05:30
temperature , &
thermalMapping , &
2013-11-27 13:34:05 +05:30
homogenization_maxNgrains , &
2014-10-09 19:38:32 +05:30
PLASTICITY_none_ID , &
2016-01-09 01:15:20 +05:30
PLASTICITY_isotropic_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_phenopowerlaw_ID , &
2017-09-30 04:06:28 +05:30
PLASTICITY_kinehardening_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_dislotwin_ID , &
2015-01-15 16:26:15 +05:30
PLASTICITY_disloucla_ID , &
2014-10-09 19:38:32 +05:30
PLASTICITY_nonlocal_ID , &
2015-05-28 22:32:23 +05:30
SOURCE_damage_isoDuctile_ID , &
SOURCE_damage_anisoBrittle_ID , &
2015-07-27 16:39:37 +05:30
SOURCE_damage_anisoDuctile_ID , &
SOURCE_thermal_externalheat_ID
2016-01-09 01:15:20 +05:30
use plastic_isotropic , only : &
plastic_isotropic_dotState
2014-12-08 21:25:30 +05:30
use plastic_phenopowerlaw , only : &
2015-01-29 19:26:09 +05:30
plastic_phenopowerlaw_dotState
2017-09-30 04:06:28 +05:30
use plastic_kinehardening , only : &
plastic_kinehardening_dotState
2014-12-08 21:25:30 +05:30
use plastic_dislotwin , only : &
2015-01-29 19:26:09 +05:30
plastic_dislotwin_dotState
2015-01-15 16:26:15 +05:30
use plastic_disloucla , only : &
2015-01-29 19:26:09 +05:30
plastic_disloucla_dotState
2014-12-08 21:25:30 +05:30
use plastic_nonlocal , only : &
2015-01-29 19:26:09 +05:30
plastic_nonlocal_dotState
2015-05-28 22:32:23 +05:30
use source_damage_isoDuctile , only : &
source_damage_isoDuctile_dotState
use source_damage_anisoBrittle , only : &
source_damage_anisoBrittle_dotState
use source_damage_anisoDuctile , only : &
source_damage_anisoDuctile_dotState
2015-07-27 16:39:37 +05:30
use source_thermal_externalheat , only : &
source_thermal_externalheat_dotState
2014-06-14 02:23:17 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) :: &
subdt !< timestep
real ( pReal ) , intent ( in ) , dimension ( homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: &
2014-03-13 12:13:49 +05:30
subfracArray !< subfraction of timestep
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: &
2014-03-13 12:13:49 +05:30
FeArray , & !< elastic deformation gradient
FpArray !< plastic deformation gradient
2018-08-25 19:29:34 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
Fi !< intermediate deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 6 ) :: &
2018-08-28 18:32:58 +05:30
S6 !< 2nd Piola Kirchhoff stress (vector notation)
2018-08-25 19:29:34 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
2018-09-17 01:46:06 +05:30
Mp
2015-05-28 22:32:23 +05:30
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
ho , & !< homogenization
tme , & !< thermal member position
2019-01-24 11:26:43 +05:30
s , & !< counter in source loop
2018-09-17 01:46:06 +05:30
instance , of
2015-10-14 00:22:01 +05:30
2018-10-04 10:09:03 +05:30
ho = material_homogenizationAt ( el )
2016-01-16 22:57:19 +05:30
tme = thermalMapping ( ho ) % p ( ip , el )
2019-01-27 15:18:27 +05:30
Mp = math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , Fi ) , math_6toSym33 ( S6 ) )
2018-08-25 19:29:34 +05:30
2016-01-16 22:57:19 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_ISOTROPIC_ID ) plasticityType
2018-12-31 16:05:01 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2018-12-30 19:44:43 +05:30
call plastic_isotropic_dotState ( Mp , instance , of )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_PHENOPOWERLAW_ID ) plasticityType
2018-09-17 01:46:06 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2018-09-17 02:27:50 +05:30
call plastic_phenopowerlaw_dotState ( Mp , instance , of )
2018-08-25 18:12:44 +05:30
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_KINEHARDENING_ID ) plasticityType
2018-12-13 14:04:40 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_kinehardening_dotState ( Mp , instance , of )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2018-10-02 01:44:54 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_dislotwin_dotState ( Mp , temperature ( ho ) % p ( tme ) , instance , of )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOUCLA_ID ) plasticityType
2018-11-30 14:34:41 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_disloucla_dotState ( Mp , temperature ( ho ) % p ( tme ) , instance , of )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2019-01-27 15:18:27 +05:30
call plastic_nonlocal_dotState ( math_sym33to6 ( Mp ) , FeArray , FpArray , temperature ( ho ) % p ( tme ) , &
2015-01-29 19:26:09 +05:30
subdt , subfracArray , ip , el )
2016-01-16 22:57:19 +05:30
end select plasticityType
2015-10-14 00:22:01 +05:30
2016-01-16 22:57:19 +05:30
SourceLoop : do s = 1_pInt , phase_Nsources ( material_phase ( ipc , ip , el ) )
2018-08-25 18:12:44 +05:30
sourceType : select case ( phase_source ( s , material_phase ( ipc , ip , el ) ) )
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoBrittle_ID ) sourceType
2018-08-28 18:32:58 +05:30
call source_damage_anisoBrittle_dotState ( S6 , ipc , ip , el ) !< correct stress?
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_isoDuctile_ID ) sourceType
2015-07-27 16:39:37 +05:30
call source_damage_isoDuctile_dotState ( ipc , ip , el )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoDuctile_ID ) sourceType
2015-07-27 16:39:37 +05:30
call source_damage_anisoDuctile_dotState ( ipc , ip , el )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
case ( SOURCE_thermal_externalheat_ID ) sourceType
2015-07-27 16:39:37 +05:30
call source_thermal_externalheat_dotState ( ipc , ip , el )
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
end select sourceType
2018-08-25 18:12:44 +05:30
2016-01-16 22:57:19 +05:30
enddo SourceLoop
2014-09-23 16:08:20 +05:30
2012-10-02 18:23:25 +05:30
end subroutine constitutive_collectDotState
2009-03-06 15:32:36 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2016-01-16 22:57:19 +05:30
!> @brief for constitutive models having an instantaneous change of state
2014-06-17 12:24:49 +05:30
!> will return false if delta state is not needed/supported by the constitutive model
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
2019-01-24 11:26:43 +05:30
subroutine constitutive_collectDeltaState ( S , Fe , Fi , ipc , ip , el )
2013-01-16 15:44:57 +05:30
use prec , only : &
2014-07-02 17:57:39 +05:30
pReal , &
2013-01-16 15:44:57 +05:30
pLongInt
use debug , only : &
debug_level , &
debug_constitutive , &
debug_levelBasic
2018-08-25 19:29:34 +05:30
use math , only : &
2019-01-24 11:26:43 +05:30
math_sym33to6 , &
2018-08-25 19:29:34 +05:30
math_mul33x33
2013-01-16 15:44:57 +05:30
use material , only : &
2018-12-13 14:04:40 +05:30
phasememberAt , &
phase_plasticityInstance , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
2015-05-28 22:32:23 +05:30
phase_source , &
phase_Nsources , &
2013-11-27 13:34:05 +05:30
material_phase , &
2017-09-30 04:06:28 +05:30
PLASTICITY_KINEHARDENING_ID , &
2015-05-28 22:32:23 +05:30
PLASTICITY_NONLOCAL_ID , &
2018-12-22 13:07:58 +05:30
SOURCE_damage_isoBrittle_ID
2017-09-30 04:06:28 +05:30
use plastic_kinehardening , only : &
plastic_kinehardening_deltaState
2014-12-08 21:25:30 +05:30
use plastic_nonlocal , only : &
plastic_nonlocal_deltaState
2015-05-28 22:32:23 +05:30
use source_damage_isoBrittle , only : &
source_damage_isoBrittle_deltaState
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2015-05-28 22:32:23 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
2019-01-24 11:26:43 +05:30
S , & !< 2nd Piola Kirchhoff stress
2018-08-25 19:29:34 +05:30
Fe , & !< elastic deformation gradient
Fi !< intermediate deformation gradient
real ( pReal ) , dimension ( 3 , 3 ) :: &
2018-12-13 14:04:40 +05:30
Mp
2015-05-28 22:32:23 +05:30
integer ( pInt ) :: &
2019-01-24 11:26:43 +05:30
i , &
2018-12-13 14:04:40 +05:30
instance , of
2013-01-16 15:44:57 +05:30
2019-01-24 11:26:43 +05:30
Mp = math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , Fi ) , S )
2018-08-25 19:29:34 +05:30
2017-09-30 04:06:28 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
2018-08-25 19:29:34 +05:30
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_KINEHARDENING_ID ) plasticityType
2018-12-13 14:04:40 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
call plastic_kinehardening_deltaState ( Mp , instance , of )
2018-08-25 19:29:34 +05:30
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2019-01-24 11:26:43 +05:30
call plastic_nonlocal_deltaState ( math_sym33to6 ( Mp ) , ip , el )
2018-08-25 19:29:34 +05:30
2017-09-30 04:06:28 +05:30
end select plasticityType
2012-05-16 20:13:26 +05:30
2019-01-24 11:26:43 +05:30
sourceLoop : do i = 1_pInt , phase_Nsources ( material_phase ( ipc , ip , el ) )
2018-08-25 19:29:34 +05:30
2019-01-24 11:26:43 +05:30
sourceType : select case ( phase_source ( i , material_phase ( ipc , ip , el ) ) )
2018-08-25 19:29:34 +05:30
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_isoBrittle_ID ) sourceType
2015-05-28 22:32:23 +05:30
call source_damage_isoBrittle_deltaState ( constitutive_homogenizedC ( ipc , ip , el ) , Fe , &
2018-08-25 19:29:34 +05:30
ipc , ip , el )
2016-01-16 22:57:19 +05:30
end select sourceType
2018-08-25 19:29:34 +05:30
2016-01-16 22:57:19 +05:30
enddo SourceLoop
2015-05-28 22:32:23 +05:30
2015-06-01 21:32:27 +05:30
end subroutine constitutive_collectDeltaState
2014-09-10 23:56:12 +05:30
2014-09-22 23:45:19 +05:30
2013-01-16 15:44:57 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns array of constitutive results
!--------------------------------------------------------------------------------------------------
2018-09-17 02:27:50 +05:30
function constitutive_postResults ( S6 , Fi , FeArray , ipc , ip , el )
2014-07-02 17:57:39 +05:30
use prec , only : &
2015-10-14 00:22:01 +05:30
pReal
2018-09-17 02:27:50 +05:30
use math , only : &
2019-01-27 15:18:27 +05:30
math_6toSym33 , &
2018-09-17 02:27:50 +05:30
math_mul33x33
2013-01-16 15:44:57 +05:30
use mesh , only : &
mesh_NcpElems , &
mesh_maxNips
use material , only : &
2018-09-17 02:27:50 +05:30
phasememberAt , &
phase_plasticityInstance , &
2014-05-22 20:54:12 +05:30
plasticState , &
2015-05-28 22:32:23 +05:30
sourceState , &
2013-01-16 15:44:57 +05:30
phase_plasticity , &
2015-05-28 22:32:23 +05:30
phase_source , &
phase_Nsources , &
2013-01-16 15:44:57 +05:30
material_phase , &
2018-10-04 10:09:03 +05:30
material_homogenizationAt , &
2015-05-28 22:32:23 +05:30
temperature , &
thermalMapping , &
2013-11-27 13:34:05 +05:30
homogenization_maxNgrains , &
PLASTICITY_NONE_ID , &
2016-01-09 01:15:20 +05:30
PLASTICITY_ISOTROPIC_ID , &
2013-11-27 13:34:05 +05:30
PLASTICITY_PHENOPOWERLAW_ID , &
2017-09-30 04:06:28 +05:30
PLASTICITY_KINEHARDENING_ID , &
2013-11-27 13:34:05 +05:30
PLASTICITY_DISLOTWIN_ID , &
2015-01-15 16:26:15 +05:30
PLASTICITY_DISLOUCLA_ID , &
2014-09-23 16:08:20 +05:30
PLASTICITY_NONLOCAL_ID , &
2015-05-28 22:32:23 +05:30
SOURCE_damage_isoBrittle_ID , &
SOURCE_damage_isoDuctile_ID , &
SOURCE_damage_anisoBrittle_ID , &
SOURCE_damage_anisoDuctile_ID
2016-01-09 01:15:20 +05:30
use plastic_isotropic , only : &
plastic_isotropic_postResults
2014-12-08 21:25:30 +05:30
use plastic_phenopowerlaw , only : &
plastic_phenopowerlaw_postResults
2017-09-30 04:06:28 +05:30
use plastic_kinehardening , only : &
plastic_kinehardening_postResults
2014-12-08 21:25:30 +05:30
use plastic_dislotwin , only : &
plastic_dislotwin_postResults
2015-01-15 16:26:15 +05:30
use plastic_disloucla , only : &
plastic_disloucla_postResults
2014-12-08 21:25:30 +05:30
use plastic_nonlocal , only : &
plastic_nonlocal_postResults
2015-05-28 22:32:23 +05:30
use source_damage_isoBrittle , only : &
source_damage_isoBrittle_postResults
use source_damage_isoDuctile , only : &
source_damage_isoDuctile_postResults
use source_damage_anisoBrittle , only : &
source_damage_anisoBrittle_postResults
use source_damage_anisoDuctile , only : &
source_damage_anisoDuctile_postResults
2014-09-23 16:08:20 +05:30
2013-01-16 15:44:57 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
2016-01-16 22:57:19 +05:30
ipc , & !< component-ID of integration point
ip , & !< integration point
el !< element
2014-09-23 16:08:20 +05:30
real ( pReal ) , dimension ( plasticState ( material_phase ( ipc , ip , el ) ) % sizePostResults + &
2015-10-14 00:22:01 +05:30
sum ( sourceState ( material_phase ( ipc , ip , el ) ) % p ( : ) % sizePostResults ) ) :: &
2014-09-23 16:08:20 +05:30
constitutive_postResults
2018-09-17 02:27:50 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
Fi !< intermediate deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , homogenization_maxNgrains , mesh_maxNips , mesh_NcpElems ) :: &
2014-03-13 12:13:49 +05:30
FeArray !< elastic deformation gradient
2013-01-16 15:44:57 +05:30
real ( pReal ) , intent ( in ) , dimension ( 6 ) :: &
2018-08-28 18:32:58 +05:30
S6 !< 2nd Piola Kirchhoff stress (vector notation)
2018-09-17 02:27:50 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
Mp !< Mandel stress
2014-10-28 08:12:25 +05:30
integer ( pInt ) :: &
2016-01-16 22:57:19 +05:30
startPos , endPos
integer ( pInt ) :: &
ho , & !< homogenization
tme , & !< thermal member position
2018-09-17 02:27:50 +05:30
s , of , instance !< counter in source loop
2015-10-14 00:22:01 +05:30
2013-01-16 15:44:57 +05:30
constitutive_postResults = 0.0_pReal
2015-10-14 00:22:01 +05:30
2019-01-27 15:18:27 +05:30
Mp = math_mul33x33 ( math_mul33x33 ( transpose ( Fi ) , Fi ) , math_6toSym33 ( S6 ) )
2018-09-17 02:27:50 +05:30
2018-10-04 10:09:03 +05:30
ho = material_homogenizationAt ( el )
2016-01-16 22:57:19 +05:30
tme = thermalMapping ( ho ) % p ( ip , el )
2015-05-28 22:32:23 +05:30
2014-09-23 16:08:20 +05:30
startPos = 1_pInt
endPos = plasticState ( material_phase ( ipc , ip , el ) ) % sizePostResults
2016-01-16 22:57:19 +05:30
plasticityType : select case ( phase_plasticity ( material_phase ( ipc , ip , el ) ) )
case ( PLASTICITY_ISOTROPIC_ID ) plasticityType
2018-12-30 19:44:43 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2018-09-17 02:27:50 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-12-30 19:44:43 +05:30
plastic_isotropic_postResults ( Mp , instance , of )
2018-10-02 01:44:54 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_PHENOPOWERLAW_ID ) plasticityType
2018-09-17 02:27:50 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2014-09-23 16:08:20 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-09-17 02:27:50 +05:30
plastic_phenopowerlaw_postResults ( Mp , instance , of )
2018-10-02 01:44:54 +05:30
2017-09-30 04:06:28 +05:30
case ( PLASTICITY_KINEHARDENING_ID ) plasticityType
2018-12-13 14:04:40 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2017-09-30 04:06:28 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-12-13 14:04:40 +05:30
plastic_kinehardening_postResults ( Mp , instance , of )
2018-10-02 01:44:54 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOTWIN_ID ) plasticityType
2018-10-02 01:44:54 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2014-09-23 16:08:20 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-10-02 01:44:54 +05:30
plastic_dislotwin_postResults ( Mp , temperature ( ho ) % p ( tme ) , instance , of )
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_DISLOUCLA_ID ) plasticityType
2018-11-30 14:34:41 +05:30
of = phasememberAt ( ipc , ip , el )
instance = phase_plasticityInstance ( material_phase ( ipc , ip , el ) )
2015-01-15 16:26:15 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-11-30 14:34:41 +05:30
plastic_disloucla_postResults ( Mp , temperature ( ho ) % p ( tme ) , instance , of )
2018-10-02 01:44:54 +05:30
2016-01-16 22:57:19 +05:30
case ( PLASTICITY_NONLOCAL_ID ) plasticityType
2014-09-23 16:08:20 +05:30
constitutive_postResults ( startPos : endPos ) = &
2018-08-28 18:32:58 +05:30
plastic_nonlocal_postResults ( S6 , FeArray , ip , el )
2016-01-16 22:57:19 +05:30
end select plasticityType
2014-09-23 16:08:20 +05:30
2016-01-16 22:57:19 +05:30
SourceLoop : do s = 1_pInt , phase_Nsources ( material_phase ( ipc , ip , el ) )
2015-05-28 22:32:23 +05:30
startPos = endPos + 1_pInt
2016-01-16 22:57:19 +05:30
endPos = endPos + sourceState ( material_phase ( ipc , ip , el ) ) % p ( s ) % sizePostResults
sourceType : select case ( phase_source ( s , material_phase ( ipc , ip , el ) ) )
case ( SOURCE_damage_isoBrittle_ID ) sourceType
2015-05-28 22:32:23 +05:30
constitutive_postResults ( startPos : endPos ) = source_damage_isoBrittle_postResults ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_isoDuctile_ID ) sourceType
2015-05-28 22:32:23 +05:30
constitutive_postResults ( startPos : endPos ) = source_damage_isoDuctile_postResults ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoBrittle_ID ) sourceType
2015-05-28 22:32:23 +05:30
constitutive_postResults ( startPos : endPos ) = source_damage_anisoBrittle_postResults ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
case ( SOURCE_damage_anisoDuctile_ID ) sourceType
2015-05-28 22:32:23 +05:30
constitutive_postResults ( startPos : endPos ) = source_damage_anisoDuctile_postResults ( ipc , ip , el )
2016-01-16 22:57:19 +05:30
end select sourceType
enddo SourceLoop
2014-10-11 02:25:09 +05:30
2012-10-02 18:23:25 +05:30
end function constitutive_postResults
2009-03-06 15:32:36 +05:30
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
2018-12-12 11:10:57 +05:30
!> @brief writes constitutive results to HDF5 output file
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
subroutine constitutive_results ( )
2018-12-12 11:10:57 +05:30
use material , only : &
PLASTICITY_ISOTROPIC_ID , &
PLASTICITY_PHENOPOWERLAW_ID , &
PLASTICITY_KINEHARDENING_ID , &
PLASTICITY_DISLOTWIN_ID , &
PLASTICITY_DISLOUCLA_ID , &
PLASTICITY_NONLOCAL_ID
#if defined(PETSc) || defined(DAMASKHDF5)
2018-12-12 12:15:20 +05:30
use results
use HDF5_utilities
use config , only : &
config_name_phase = > phase_name ! anticipate logical name
2018-12-13 15:21:35 +05:30
2018-12-12 12:15:20 +05:30
use material , only : &
2018-12-13 15:21:35 +05:30
phase_plasticityInstance , &
2018-12-12 12:15:20 +05:30
material_phase_plasticity_type = > phase_plasticity
2018-12-12 11:10:57 +05:30
use plastic_phenopowerlaw , only : &
plastic_phenopowerlaw_results
2018-12-12 12:15:20 +05:30
implicit none
integer ( pInt ) :: p
call HDF5_closeGroup ( results_addGroup ( 'current/phase' ) )
do p = 1 , size ( config_name_phase )
call HDF5_closeGroup ( results_addGroup ( 'current/phase/' / / trim ( config_name_phase ( p ) ) ) )
if ( material_phase_plasticity_type ( p ) == PLASTICITY_PHENOPOWERLAW_ID ) then
2018-12-13 15:21:35 +05:30
call plastic_phenopowerlaw_results ( phase_plasticityInstance ( p ) , 'current/phase/' / / trim ( config_name_phase ( p ) ) )
2018-12-12 12:15:20 +05:30
endif
enddo
2018-12-12 11:10:57 +05:30
#endif
2018-12-05 04:25:39 +05:30
end subroutine constitutive_results
2013-02-11 16:13:45 +05:30
end module constitutive