2012-10-02 18:23:25 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-03-06 20:11:15 +05:30
|
|
|
! $Id$
|
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
|
2013-01-16 15:44:57 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
private
|
2014-05-27 20:16:03 +05:30
|
|
|
integer(pInt), public, protected :: &
|
|
|
|
constitutive_maxSizePostResults, &
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_maxSizeDotState, &
|
|
|
|
constitutive_damage_maxSizePostResults, &
|
|
|
|
constitutive_damage_maxSizeDotState, &
|
|
|
|
constitutive_thermal_maxSizePostResults, &
|
2014-10-11 02:25:09 +05:30
|
|
|
constitutive_thermal_maxSizeDotState, &
|
|
|
|
constitutive_vacancy_maxSizePostResults, &
|
|
|
|
constitutive_vacancy_maxSizeDotState
|
2014-06-26 19:23:12 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
public :: &
|
|
|
|
constitutive_init, &
|
|
|
|
constitutive_homogenizedC, &
|
2014-10-27 21:03:35 +05:30
|
|
|
constitutive_damagedC, &
|
2013-01-16 15:44:57 +05:30
|
|
|
constitutive_microstructure, &
|
|
|
|
constitutive_LpAndItsTangent, &
|
2014-11-01 00:33:08 +05:30
|
|
|
constitutive_LiAndItsTangent, &
|
|
|
|
constitutive_getFi, &
|
2014-11-03 16:13:36 +05:30
|
|
|
constitutive_putFi, &
|
2014-11-01 00:33:08 +05:30
|
|
|
constitutive_getFi0, &
|
|
|
|
constitutive_getPartionedFi0, &
|
2013-01-16 15:44:57 +05:30
|
|
|
constitutive_TandItsTangent, &
|
|
|
|
constitutive_collectDotState, &
|
|
|
|
constitutive_collectDeltaState, &
|
2014-09-05 22:01:27 +05:30
|
|
|
constitutive_getLocalDamage, &
|
2014-09-22 23:45:19 +05:30
|
|
|
constitutive_putLocalDamage, &
|
2014-09-26 23:37:48 +05:30
|
|
|
constitutive_getDamage, &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_getSlipDamage, &
|
2014-10-11 15:39:36 +05:30
|
|
|
constitutive_getDamageDiffusion33, &
|
2014-09-26 23:37:48 +05:30
|
|
|
constitutive_getAdiabaticTemperature, &
|
|
|
|
constitutive_putAdiabaticTemperature, &
|
|
|
|
constitutive_getTemperature, &
|
2014-10-11 02:25:09 +05:30
|
|
|
constitutive_getLocalVacancyConcentration, &
|
|
|
|
constitutive_putLocalVacancyConcentration, &
|
|
|
|
constitutive_getVacancyConcentration, &
|
2014-10-11 16:09:44 +05:30
|
|
|
constitutive_getVacancyDiffusion33, &
|
2013-01-16 15:44:57 +05:30
|
|
|
constitutive_postResults
|
|
|
|
|
|
|
|
private :: &
|
2014-10-10 22:04:51 +05:30
|
|
|
constitutive_hooke_TandItsTangent, &
|
2014-10-11 15:15:30 +05:30
|
|
|
constitutive_getAccumulatedSlip, &
|
|
|
|
constitutive_getSlipRate
|
2013-01-16 15:44:57 +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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-13 23:24:27 +05:30
|
|
|
subroutine constitutive_init(temperature_init)
|
2014-03-12 13:03:51 +05:30
|
|
|
#ifdef HDF
|
|
|
|
use hdf5, only: &
|
|
|
|
HID_T
|
|
|
|
use IO, only : &
|
2014-03-12 22:21:01 +05:30
|
|
|
HDF5_mappingConstitutive
|
2014-03-12 13:03:51 +05:30
|
|
|
#endif
|
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
|
2014-07-02 17:57:39 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-01-16 15:44:57 +05:30
|
|
|
use debug, only: &
|
|
|
|
debug_level, &
|
|
|
|
debug_constitutive, &
|
|
|
|
debug_levelBasic
|
|
|
|
use numerics, only: &
|
2014-10-10 01:53:06 +05:30
|
|
|
worldrank, &
|
2013-01-16 15:44:57 +05:30
|
|
|
numerics_integrator
|
|
|
|
use IO, only: &
|
|
|
|
IO_error, &
|
|
|
|
IO_open_file, &
|
|
|
|
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
|
2013-01-16 15:44:57 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_maxNips, &
|
|
|
|
mesh_NcpElems, &
|
|
|
|
mesh_element, &
|
|
|
|
FE_Nips, &
|
|
|
|
FE_geomtype
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
material_Nphase, &
|
|
|
|
material_localFileExt, &
|
|
|
|
material_configFile, &
|
|
|
|
phase_name, &
|
|
|
|
phase_elasticity, &
|
|
|
|
phase_plasticity, &
|
|
|
|
phase_plasticityInstance, &
|
2014-09-23 16:08:20 +05:30
|
|
|
phase_damage, &
|
|
|
|
phase_damageInstance, &
|
|
|
|
phase_thermal, &
|
|
|
|
phase_thermalInstance, &
|
2014-10-11 02:25:09 +05:30
|
|
|
phase_vacancy, &
|
|
|
|
phase_vacancyInstance, &
|
2013-01-16 15:44:57 +05:30
|
|
|
phase_Noutput, &
|
|
|
|
homogenization_Ngrains, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
2014-10-09 19:38:32 +05:30
|
|
|
ELASTICITY_hooke_ID, &
|
|
|
|
PLASTICITY_none_ID, &
|
|
|
|
PLASTICITY_j2_ID, &
|
|
|
|
PLASTICITY_phenopowerlaw_ID, &
|
|
|
|
PLASTICITY_dislotwin_ID, &
|
|
|
|
PLASTICITY_dislokmc_ID, &
|
|
|
|
PLASTICITY_titanmod_ID, &
|
|
|
|
PLASTICITY_nonlocal_ID ,&
|
2013-11-27 13:34:05 +05:30
|
|
|
ELASTICITY_HOOKE_label, &
|
|
|
|
PLASTICITY_NONE_label, &
|
|
|
|
PLASTICITY_J2_label, &
|
|
|
|
PLASTICITY_PHENOPOWERLAW_label, &
|
|
|
|
PLASTICITY_DISLOTWIN_label, &
|
2014-08-08 16:34:40 +05:30
|
|
|
PLASTICITY_DISLOKMC_label, &
|
2013-11-27 13:34:05 +05:30
|
|
|
PLASTICITY_TITANMOD_label, &
|
2014-08-08 16:34:40 +05:30
|
|
|
PLASTICITY_NONLOCAL_label, &
|
2014-10-09 19:38:32 +05:30
|
|
|
LOCAL_DAMAGE_none_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
2014-10-09 19:38:32 +05:30
|
|
|
LOCAL_THERMAL_isothermal_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID, &
|
2014-10-11 02:25:09 +05:30
|
|
|
LOCAL_VACANCY_constant_ID, &
|
|
|
|
LOCAL_VACANCY_generation_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_none_LABEL, &
|
|
|
|
LOCAL_DAMAGE_isoBrittle_LABEL, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_LABEL, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_LABEL, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_LABEL, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_LABEL, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_THERMAL_isothermal_label, &
|
|
|
|
LOCAL_THERMAL_adiabatic_label, &
|
2014-10-11 02:25:09 +05:30
|
|
|
LOCAL_VACANCY_constant_label, &
|
|
|
|
LOCAL_VACANCY_generation_label, &
|
2014-05-08 20:25:19 +05:30
|
|
|
plasticState, &
|
2014-09-23 16:08:20 +05:30
|
|
|
damageState, &
|
|
|
|
thermalState, &
|
2014-10-11 02:25:09 +05:30
|
|
|
vacancyState, &
|
2014-08-08 16:34:40 +05:30
|
|
|
mappingConstitutive
|
2014-07-02 17:57:39 +05:30
|
|
|
|
2014-08-08 16:34:40 +05:30
|
|
|
|
2012-07-03 16:46:38 +05:30
|
|
|
use constitutive_none
|
2012-03-09 01:55:28 +05:30
|
|
|
use constitutive_j2
|
|
|
|
use constitutive_phenopowerlaw
|
2014-06-03 19:16:42 +05:30
|
|
|
use constitutive_dislotwin
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc
|
2012-03-09 01:55:28 +05:30
|
|
|
use constitutive_titanmod
|
2014-06-14 02:23:17 +05:30
|
|
|
use constitutive_nonlocal
|
2014-09-23 16:08:20 +05:30
|
|
|
use damage_none
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle
|
|
|
|
use damage_isoDuctile
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_anisoBrittle
|
2014-10-10 18:12:12 +05:30
|
|
|
use damage_gurson
|
2014-09-26 21:37:26 +05:30
|
|
|
use thermal_isothermal
|
2014-10-09 19:38:32 +05:30
|
|
|
use thermal_adiabatic
|
2014-10-11 02:25:09 +05:30
|
|
|
use vacancy_constant
|
|
|
|
use vacancy_generation
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
2014-10-13 23:24:27 +05:30
|
|
|
real(pReal), intent(in) :: temperature_init !< initial temperature
|
2013-12-12 22:39:59 +05:30
|
|
|
integer(pInt), parameter :: FILEUNIT = 200_pInt
|
2013-01-16 15:44:57 +05:30
|
|
|
integer(pInt) :: &
|
2014-10-13 23:24:27 +05:30
|
|
|
e, & !< maximum number of elements
|
2014-02-28 15:48:40 +05:30
|
|
|
phase, &
|
2014-07-03 18:47:29 +05:30
|
|
|
instance
|
2014-05-08 20:25:19 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
integer(pInt), dimension(:,:), pointer :: thisSize
|
2014-09-26 15:55:26 +05:30
|
|
|
integer(pInt), dimension(:) , pointer :: thisNoutput
|
2013-01-16 15:44:57 +05:30
|
|
|
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
|
2014-10-11 02:25:09 +05:30
|
|
|
logical :: knownPlasticity, knownDamage, knownThermal, knownVacancy, nonlocalConstitutionPresent
|
2013-05-24 19:13:44 +05:30
|
|
|
nonlocalConstitutionPresent = .false.
|
2013-01-16 15:44:57 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse plasticities from config file
|
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
|
2014-11-07 17:45:28 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_NONE_ID)) call constitutive_none_init
|
2014-03-09 02:20:31 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_J2_ID)) call constitutive_j2_init(FILEUNIT)
|
|
|
|
if (any(phase_plasticity == PLASTICITY_PHENOPOWERLAW_ID)) call constitutive_phenopowerlaw_init(FILEUNIT)
|
2014-06-03 19:16:42 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_DISLOTWIN_ID)) call constitutive_dislotwin_init(FILEUNIT)
|
2014-08-08 16:34:40 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_DISLOKMC_ID)) call constitutive_dislokmc_init(FILEUNIT)
|
2014-03-09 02:20:31 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_TITANMOD_ID)) call constitutive_titanmod_init(FILEUNIT)
|
2014-07-08 20:28:23 +05:30
|
|
|
if (any(phase_plasticity == PLASTICITY_NONLOCAL_ID)) then
|
|
|
|
call constitutive_nonlocal_init(FILEUNIT)
|
|
|
|
call constitutive_nonlocal_stateInit()
|
|
|
|
endif
|
2013-12-12 22:39:59 +05:30
|
|
|
close(FILEUNIT)
|
2014-09-23 16:08:20 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse damage from config file
|
|
|
|
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
|
2014-11-07 17:45:28 +05:30
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_none_ID)) call damage_none_init
|
2014-10-28 16:19:12 +05:30
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_isoBrittle_ID)) call damage_isoBrittle_init(FILEUNIT)
|
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_isoductile_ID)) call damage_isoDuctile_init(FILEUNIT)
|
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_anisoBrittle_ID)) call damage_anisoBrittle_init(FILEUNIT)
|
2014-11-05 23:11:08 +05:30
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_anisoductile_ID)) call damage_anisoDuctile_init(FILEUNIT)
|
|
|
|
if (any(phase_damage == LOCAL_DAMAGE_gurson_ID)) call damage_gurson_init(FILEUNIT)
|
2014-09-23 16:08:20 +05:30
|
|
|
close(FILEUNIT)
|
2013-01-16 15:44:57 +05:30
|
|
|
|
2014-09-23 16:08:20 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse thermal from config file
|
|
|
|
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
|
2014-11-07 17:45:28 +05:30
|
|
|
if (any(phase_thermal == LOCAL_THERMAL_isothermal_ID)) call thermal_isothermal_init(temperature_init)
|
2014-10-13 23:24:27 +05:30
|
|
|
if (any(phase_thermal == LOCAL_THERMAL_adiabatic_ID)) call thermal_adiabatic_init(FILEUNIT,temperature_init)
|
2014-10-11 02:25:09 +05:30
|
|
|
close(FILEUNIT)
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse vacancy model from config file
|
|
|
|
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
|
2014-11-07 17:45:28 +05:30
|
|
|
if (any(phase_vacancy == LOCAL_VACANCY_constant_ID)) call vacancy_constant_init
|
2014-10-11 02:25:09 +05:30
|
|
|
if (any(phase_vacancy == LOCAL_VACANCY_generation_ID)) call vacancy_generation_init(FILEUNIT)
|
2014-09-23 16:08:20 +05:30
|
|
|
close(FILEUNIT)
|
|
|
|
|
2014-10-10 22:15:14 +05:30
|
|
|
mainProcess: if (worldrank == 0) then
|
|
|
|
write(6,'(/,a)') ' <<<+- constitutive init -+>>>'
|
|
|
|
write(6,'(a)') ' $Id$'
|
|
|
|
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
|
2012-10-09 18:04:57 +05:30
|
|
|
#include "compilation_info.f90"
|
2014-10-10 22:15:14 +05:30
|
|
|
endif mainProcess
|
2013-01-16 15:44:57 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! write description file for constitutive phase output
|
2013-12-12 22:39:59 +05:30
|
|
|
call IO_write_jobFile(FILEUNIT,'outputConstitutive')
|
2014-02-28 15:48:40 +05:30
|
|
|
do phase = 1_pInt,material_Nphase
|
|
|
|
instance = phase_plasticityInstance(phase) ! which instance of a plasticity is present phase
|
2013-01-16 15:44:57 +05:30
|
|
|
knownPlasticity = .true. ! assume valid
|
2014-02-28 15:48:40 +05:30
|
|
|
select case(phase_plasticity(phase)) ! split per constititution
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONE_ID)
|
|
|
|
outputName = PLASTICITY_NONE_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => null()
|
2014-03-13 12:13:49 +05:30
|
|
|
thisOutput => null() ! constitutive_none_output
|
|
|
|
thisSize => null() ! constitutive_none_sizePostResult
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_J2_ID)
|
|
|
|
outputName = PLASTICITY_J2_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_j2_Noutput
|
2013-01-16 15:44:57 +05:30
|
|
|
thisOutput => constitutive_j2_output
|
|
|
|
thisSize => constitutive_j2_sizePostResult
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
2013-12-18 12:58:01 +05:30
|
|
|
outputName = PLASTICITY_PHENOPOWERLAW_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_phenopowerlaw_Noutput
|
2013-01-16 15:44:57 +05:30
|
|
|
thisOutput => constitutive_phenopowerlaw_output
|
|
|
|
thisSize => constitutive_phenopowerlaw_sizePostResult
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
|
|
|
outputName = PLASTICITY_DISLOTWIN_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_dislotwin_Noutput
|
2013-01-16 15:44:57 +05:30
|
|
|
thisOutput => constitutive_dislotwin_output
|
|
|
|
thisSize => constitutive_dislotwin_sizePostResult
|
2014-08-08 16:34:40 +05:30
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
|
|
|
outputName = PLASTICITY_DISLOKMC_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_dislokmc_Noutput
|
2014-08-08 16:34:40 +05:30
|
|
|
thisOutput => constitutive_dislokmc_output
|
|
|
|
thisSize => constitutive_dislokmc_sizePostResult
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
|
|
|
outputName = PLASTICITY_TITANMOD_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_titanmod_Noutput
|
2013-11-27 13:34:05 +05:30
|
|
|
thisOutput => constitutive_titanmod_output
|
|
|
|
thisSize => constitutive_titanmod_sizePostResult
|
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
|
|
|
outputName = PLASTICITY_NONLOCAL_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => constitutive_nonlocal_Noutput
|
2013-01-16 15:44:57 +05:30
|
|
|
thisOutput => constitutive_nonlocal_output
|
|
|
|
thisSize => constitutive_nonlocal_sizePostResult
|
|
|
|
case default
|
|
|
|
knownPlasticity = .false.
|
|
|
|
end select
|
2014-02-28 15:48:40 +05:30
|
|
|
write(FILEUNIT,'(/,a,/)') '['//trim(phase_name(phase))//']'
|
2013-01-16 15:44:57 +05:30
|
|
|
if (knownPlasticity) then
|
2013-12-12 22:39:59 +05:30
|
|
|
write(FILEUNIT,'(a)') '(plasticity)'//char(9)//trim(outputName)
|
2014-06-25 04:52:52 +05:30
|
|
|
if (phase_plasticity(phase) /= PLASTICITY_NONE_ID) then
|
2014-09-26 15:55:26 +05:30
|
|
|
do e = 1_pInt,thisNoutput(instance)
|
2014-06-25 04:52:52 +05:30
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,instance))//char(9),thisSize(e,instance)
|
|
|
|
enddo
|
|
|
|
endif
|
2013-01-16 15:44:57 +05:30
|
|
|
endif
|
2014-09-23 16:08:20 +05:30
|
|
|
#ifdef multiphysicsOut
|
|
|
|
instance = phase_damageInstance(phase) ! which instance of a plasticity is present phase
|
|
|
|
knownDamage = .true.
|
|
|
|
select case(phase_damage(phase)) ! split per constititution
|
|
|
|
case (LOCAL_DAMAGE_none_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_NONE_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => null()
|
2014-09-23 16:08:20 +05:30
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_isoBrittle_LABEL
|
|
|
|
thisNoutput => damage_isoBrittle_Noutput
|
|
|
|
thisOutput => damage_isoBrittle_output
|
|
|
|
thisSize => damage_isoBrittle_sizePostResult
|
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_isoDuctile_LABEL
|
|
|
|
thisNoutput => damage_isoDuctile_Noutput
|
|
|
|
thisOutput => damage_isoDuctile_output
|
|
|
|
thisSize => damage_isoDuctile_sizePostResult
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_anisoBrittle_label
|
|
|
|
thisNoutput => damage_anisoBrittle_Noutput
|
|
|
|
thisOutput => damage_anisoBrittle_output
|
|
|
|
thisSize => damage_anisoBrittle_sizePostResult
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_anisoDuctile_LABEL
|
|
|
|
thisNoutput => damage_anisoDuctile_Noutput
|
|
|
|
thisOutput => damage_anisoDuctile_output
|
|
|
|
thisSize => damage_anisoDuctile_sizePostResult
|
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
outputName = LOCAL_DAMAGE_gurson_label
|
|
|
|
thisNoutput => damage_gurson_Noutput
|
|
|
|
thisOutput => damage_gurson_output
|
|
|
|
thisSize => damage_gurson_sizePostResult
|
2014-09-23 16:08:20 +05:30
|
|
|
case default
|
|
|
|
knownDamage = .false.
|
|
|
|
end select
|
|
|
|
if (knownDamage) then
|
|
|
|
write(FILEUNIT,'(a)') '(damage)'//char(9)//trim(outputName)
|
|
|
|
if (phase_damage(phase) /= LOCAL_DAMAGE_none_ID) then
|
2014-09-26 15:55:26 +05:30
|
|
|
do e = 1_pInt,thisNoutput(instance)
|
2014-09-23 16:08:20 +05:30
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,instance))//char(9),thisSize(e,instance)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
instance = phase_thermalInstance(phase) ! which instance is present phase
|
|
|
|
knownThermal = .true.
|
|
|
|
select case(phase_thermal(phase)) ! split per constititution
|
2014-10-11 02:25:09 +05:30
|
|
|
case (LOCAL_THERMAL_isothermal_ID)
|
2014-09-26 21:37:26 +05:30
|
|
|
outputName = LOCAL_THERMAL_ISOTHERMAL_label
|
2014-09-26 15:55:26 +05:30
|
|
|
thisNoutput => null()
|
2014-09-23 16:08:20 +05:30
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
2014-10-09 19:38:32 +05:30
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
outputName = LOCAL_THERMAL_ADIABATIC_label
|
|
|
|
thisNoutput => thermal_adiabatic_Noutput
|
|
|
|
thisOutput => thermal_adiabatic_output
|
|
|
|
thisSize => thermal_adiabatic_sizePostResult
|
2014-09-23 16:08:20 +05:30
|
|
|
case default
|
|
|
|
knownThermal = .false.
|
|
|
|
end select
|
|
|
|
if (knownThermal) then
|
|
|
|
write(FILEUNIT,'(a)') '(thermal)'//char(9)//trim(outputName)
|
2014-10-09 19:38:32 +05:30
|
|
|
if (phase_thermal(phase) /= LOCAL_THERMAL_isothermal_ID) then
|
2014-09-26 15:55:26 +05:30
|
|
|
do e = 1_pInt,thisNoutput(instance)
|
2014-09-23 16:08:20 +05:30
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,instance))//char(9),thisSize(e,instance)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
endif
|
2014-10-11 02:25:09 +05:30
|
|
|
instance = phase_vacancyInstance(phase) ! which instance is present phase
|
|
|
|
knownVacancy = .true.
|
|
|
|
select case(phase_vacancy(phase)) ! split per constititution
|
|
|
|
case (LOCAL_VACANCY_constant_ID)
|
|
|
|
outputName = LOCAL_VACANCY_constant_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
|
|
|
outputName = LOCAL_VACANCY_generation_label
|
|
|
|
thisNoutput => vacancy_generation_Noutput
|
|
|
|
thisOutput => vacancy_generation_output
|
|
|
|
thisSize => vacancy_generation_sizePostResult
|
|
|
|
case default
|
|
|
|
knownVacancy = .false.
|
|
|
|
end select
|
|
|
|
if (knownVacancy) then
|
|
|
|
write(FILEUNIT,'(a)') '(vacancy)'//char(9)//trim(outputName)
|
2014-10-13 23:24:27 +05:30
|
|
|
if (phase_vacancy(phase) /= LOCAL_VACANCY_constant_ID) then
|
2014-10-11 02:25:09 +05:30
|
|
|
do e = 1_pInt,thisNoutput(instance)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,instance))//char(9),thisSize(e,instance)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
endif
|
2014-09-23 16:08:20 +05:30
|
|
|
#endif
|
2013-01-16 15:44:57 +05:30
|
|
|
enddo
|
2013-12-12 22:39:59 +05:30
|
|
|
close(FILEUNIT)
|
2013-01-16 15:44:57 +05:30
|
|
|
|
2014-07-03 18:47:29 +05:30
|
|
|
constitutive_maxSizeDotState = 0_pInt
|
|
|
|
constitutive_maxSizePostResults = 0_pInt
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_damage_maxSizePostResults = 0_pInt
|
|
|
|
constitutive_damage_maxSizeDotState = 0_pInt
|
|
|
|
constitutive_thermal_maxSizePostResults = 0_pInt
|
|
|
|
constitutive_thermal_maxSizeDotState = 0_pInt
|
2014-10-11 02:25:09 +05:30
|
|
|
constitutive_vacancy_maxSizePostResults = 0_pInt
|
|
|
|
constitutive_vacancy_maxSizeDotState = 0_pInt
|
2014-07-03 18:47:29 +05:30
|
|
|
|
|
|
|
PhaseLoop2:do phase = 1_pInt,material_Nphase
|
|
|
|
plasticState(phase)%partionedState0 = plasticState(phase)%State0
|
|
|
|
plasticState(phase)%State = plasticState(phase)%State0
|
|
|
|
constitutive_maxSizeDotState = max(constitutive_maxSizeDotState, plasticState(phase)%sizeDotState)
|
|
|
|
constitutive_maxSizePostResults = max(constitutive_maxSizePostResults, plasticState(phase)%sizePostResults)
|
2014-09-23 16:08:20 +05:30
|
|
|
damageState(phase)%partionedState0 = damageState(phase)%State0
|
|
|
|
damageState(phase)%State = damageState(phase)%State0
|
|
|
|
constitutive_damage_maxSizeDotState = max(constitutive_damage_maxSizeDotState, damageState(phase)%sizeDotState)
|
|
|
|
constitutive_damage_maxSizePostResults = max(constitutive_damage_maxSizePostResults, damageState(phase)%sizePostResults)
|
|
|
|
thermalState(phase)%partionedState0 = thermalState(phase)%State0
|
|
|
|
thermalState(phase)%State = thermalState(phase)%State0
|
|
|
|
constitutive_thermal_maxSizeDotState = max(constitutive_thermal_maxSizeDotState, thermalState(phase)%sizeDotState)
|
|
|
|
constitutive_thermal_maxSizePostResults = max(constitutive_thermal_maxSizePostResults, thermalState(phase)%sizePostResults)
|
2014-10-11 02:25:09 +05:30
|
|
|
vacancyState(phase)%partionedState0 = vacancyState(phase)%State0
|
|
|
|
vacancyState(phase)%State = vacancyState(phase)%State0
|
|
|
|
constitutive_vacancy_maxSizeDotState = max(constitutive_vacancy_maxSizeDotState, vacancyState(phase)%sizeDotState)
|
|
|
|
constitutive_vacancy_maxSizePostResults = max(constitutive_vacancy_maxSizePostResults, vacancyState(phase)%sizePostResults)
|
2014-07-03 18:47:29 +05:30
|
|
|
enddo PhaseLoop2
|
2014-06-24 14:54:59 +05:30
|
|
|
|
2014-03-12 13:03:51 +05:30
|
|
|
#ifdef HDF
|
2014-03-26 14:11:45 +05:30
|
|
|
call HDF5_mappingConstitutive(mappingConstitutive)
|
2014-04-15 15:39:20 +05:30
|
|
|
do phase = 1_pInt,material_Nphase
|
|
|
|
instance = phase_plasticityInstance(phase) ! which instance of a plasticity is present phase
|
|
|
|
select case(phase_plasticity(phase)) ! split per constititution
|
|
|
|
case (PLASTICITY_NONE_ID)
|
|
|
|
case (PLASTICITY_J2_ID)
|
|
|
|
end select
|
|
|
|
enddo
|
2014-03-12 13:03:51 +05:30
|
|
|
#endif
|
2014-05-08 20:25:19 +05:30
|
|
|
|
2014-07-02 17:57:39 +05:30
|
|
|
#ifdef TODO
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! report
|
|
|
|
constitutive_maxSizeState = maxval(constitutive_sizeState)
|
|
|
|
constitutive_maxSizeDotState = maxval(constitutive_sizeDotState)
|
|
|
|
|
|
|
|
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0_pInt) then
|
2013-11-27 13:34:05 +05:30
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_state0: ', shape(constitutive_state0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_partionedState0: ', shape(constitutive_partionedState0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_subState0: ', shape(constitutive_subState0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_state: ', shape(constitutive_state)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_aTolState: ', shape(constitutive_aTolState)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_dotState: ', shape(constitutive_dotState)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_deltaState: ', shape(constitutive_deltaState)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_sizeState: ', shape(constitutive_sizeState)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'constitutive_sizeDotState: ', shape(constitutive_sizeDotState)
|
|
|
|
write(6,'(a32,1x,7(i8,1x),/)') 'constitutive_sizePostResults: ', shape(constitutive_sizePostResults)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'maxSizeState: ', constitutive_maxSizeState
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'maxSizeDotState: ', constitutive_maxSizeDotState
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'maxSizePostResults: ', constitutive_maxSizePostResults
|
2013-01-16 15:44:57 +05:30
|
|
|
endif
|
|
|
|
flush(6)
|
2014-07-02 17:57:39 +05:30
|
|
|
#endif
|
|
|
|
|
|
|
|
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-07-02 17:57:39 +05:30
|
|
|
function constitutive_homogenizedC(ipc,ip,el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-01-16 15:44:57 +05:30
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
2013-11-27 13:34:05 +05:30
|
|
|
material_phase, &
|
|
|
|
PLASTICITY_TITANMOD_ID, &
|
2014-08-08 16:34:40 +05:30
|
|
|
PLASTICITY_DISLOTWIN_ID, &
|
|
|
|
PLASTICITY_DISLOKMC_ID, &
|
2014-06-03 19:16:42 +05:30
|
|
|
plasticState,&
|
2014-08-08 16:34:40 +05:30
|
|
|
mappingConstitutive
|
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_homogenizedC
|
2014-06-03 19:16:42 +05:30
|
|
|
use constitutive_dislotwin, only: &
|
|
|
|
constitutive_dislotwin_homogenizedC
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_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) :: &
|
2013-09-19 13:16:01 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2009-03-06 15:32:36 +05:30
|
|
|
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2014-06-03 19:16:42 +05:30
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
2014-08-08 16:34:40 +05:30
|
|
|
constitutive_homogenizedC = constitutive_dislotwin_homogenizedC(ipc,ip,el)
|
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
|
|
|
constitutive_homogenizedC = constitutive_dislokmc_homogenizedC(ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
2014-07-02 17:57:39 +05:30
|
|
|
constitutive_homogenizedC = constitutive_titanmod_homogenizedC (ipc,ip,el)
|
2014-03-09 02:20:31 +05:30
|
|
|
case default
|
2014-07-02 17:57:39 +05:30
|
|
|
constitutive_homogenizedC = lattice_C66(1:6,1:6,material_phase (ipc,ip,el))
|
2009-08-11 22:01:57 +05:30
|
|
|
|
2009-03-04 19:31:36 +05:30
|
|
|
end select
|
2009-03-06 15:32:36 +05:30
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
end function constitutive_homogenizedC
|
|
|
|
|
2014-10-27 21:03:35 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-28 06:48:10 +05:30
|
|
|
!> @brief returns the damaged elasticity matrix if relevant
|
2014-10-27 21:03:35 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_damagedC(ipc,ip,el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
2014-10-27 21:03:35 +05:30
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
2014-11-12 04:16:41 +05:30
|
|
|
damage_isoBrittle_getDamagedC66
|
2014-10-27 21:03:35 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
real(pReal), dimension(6,6) :: constitutive_damagedC
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
2014-10-28 06:48:10 +05:30
|
|
|
el
|
2014-10-27 21:03:35 +05:30
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
2014-11-12 04:16:41 +05:30
|
|
|
constitutive_damagedC = damage_isoBrittle_getDamagedC66(constitutive_homogenizedC(ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2014-10-28 06:48:10 +05:30
|
|
|
case default
|
|
|
|
constitutive_damagedC = constitutive_homogenizedC(ipc,ip,el)
|
2014-10-27 21:03:35 +05:30
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_damagedC
|
2009-03-06 15:32:36 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief calls microstructure function of the different constitutive models
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-10 17:58:57 +05:30
|
|
|
subroutine constitutive_microstructure(Tstar_v, Fe, Fp, ipc, ip, el)
|
2014-07-02 17:57:39 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-01-16 15:44:57 +05:30
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
2014-09-23 16:08:20 +05:30
|
|
|
phase_damage, &
|
2013-11-27 13:34:05 +05:30
|
|
|
material_phase, &
|
2014-10-10 18:12:12 +05:30
|
|
|
PLASTICITY_dislotwin_ID, &
|
|
|
|
PLASTICITY_dislokmc_ID, &
|
|
|
|
PLASTICITY_titanmod_ID, &
|
|
|
|
PLASTICITY_nonlocal_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_gurson_ID
|
2014-08-08 16:34:40 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_microstructure
|
|
|
|
use constitutive_nonlocal, only: &
|
|
|
|
constitutive_nonlocal_microstructure
|
2014-06-03 19:16:42 +05:30
|
|
|
use constitutive_dislotwin, only: &
|
|
|
|
constitutive_dislotwin_microstructure
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_microstructure
|
2014-10-10 18:12:12 +05:30
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_microstructure
|
2014-06-03 19:16:42 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2014-11-03 16:13:36 +05:30
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
2014-09-23 16:08:20 +05:30
|
|
|
Tstar_v !< 2nd Piola Kirchhoff stress tensor (Mandel)
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Fe, & !< elastic deformation gradient
|
|
|
|
Fp !< plastic deformation gradient
|
2014-10-10 22:04:51 +05:30
|
|
|
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2013-11-27 13:34:05 +05:30
|
|
|
|
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_dislotwin_microstructure(constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2014-08-08 16:34:40 +05:30
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_dislokmc_microstructure(constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_titanmod_microstructure (constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
2014-07-02 17:57:39 +05:30
|
|
|
call constitutive_nonlocal_microstructure (Fe,Fp, ip,el)
|
|
|
|
|
2014-06-03 19:16:42 +05:30
|
|
|
end select
|
2014-09-23 16:08:20 +05:30
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-10 18:12:12 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
call damage_gurson_microstructure(ipc, ip, el)
|
2014-09-23 16:08:20 +05:30
|
|
|
|
|
|
|
end select
|
2014-06-03 19:16:42 +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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-09-19 13:16:01 +05:30
|
|
|
!> @brief contains the constitutive equation for calculating the velocity gradient
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-10 17:58:57 +05:30
|
|
|
subroutine constitutive_LpAndItsTangent(Lp, dLp_dTstar, Tstar_v, ipc, ip, el)
|
2014-07-02 17:57:39 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-10-16 18:34:59 +05:30
|
|
|
use math, only: &
|
|
|
|
math_identity2nd
|
2013-01-16 15:44:57 +05:30
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
2014-10-28 08:12:25 +05:30
|
|
|
phase_plasticityInstance, &
|
2013-11-27 13:34:05 +05:30
|
|
|
material_phase, &
|
2014-05-22 20:54:12 +05:30
|
|
|
plasticState,&
|
2014-06-23 00:28:29 +05:30
|
|
|
mappingConstitutive, &
|
2013-11-27 13:34:05 +05:30
|
|
|
PLASTICITY_NONE_ID, &
|
|
|
|
PLASTICITY_J2_ID, &
|
|
|
|
PLASTICITY_PHENOPOWERLAW_ID, &
|
|
|
|
PLASTICITY_DISLOTWIN_ID, &
|
2014-08-08 16:34:40 +05:30
|
|
|
PLASTICITY_DISLOKMC_ID, &
|
2013-11-27 13:34:05 +05:30
|
|
|
PLASTICITY_TITANMOD_ID, &
|
|
|
|
PLASTICITY_NONLOCAL_ID
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_j2, only: &
|
|
|
|
constitutive_j2_LpAndItsTangent
|
|
|
|
use constitutive_phenopowerlaw, only: &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_phenopowerlaw_LpAndItsTangent, &
|
|
|
|
constitutive_phenopowerlaw_totalNslip
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_dislotwin, only: &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_dislotwin_LpAndItsTangent, &
|
|
|
|
constitutive_dislotwin_totalNslip
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc, only: &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_dislokmc_LpAndItsTangent, &
|
|
|
|
constitutive_dislokmc_totalNslip
|
2013-11-27 13:34:05 +05:30
|
|
|
use constitutive_titanmod, only: &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_titanmod_LpAndItsTangent, &
|
|
|
|
constitutive_titanmod_totalNslip
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_nonlocal, only: &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_nonlocal_LpAndItsTangent, &
|
|
|
|
totalNslip
|
2014-06-14 02:23:17 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola-Kirchhoff stress
|
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
|
|
|
Lp !< plastic velocity gradient
|
|
|
|
real(pReal), intent(out), dimension(9,9) :: &
|
|
|
|
dLp_dTstar !< derivative of Lp with respect to Tstar (4th-order tensor)
|
2014-10-28 08:12:25 +05:30
|
|
|
integer(pInt) :: &
|
|
|
|
nSlip
|
|
|
|
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2013-01-16 15:44:57 +05:30
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONE_ID)
|
2013-10-16 18:34:59 +05:30
|
|
|
Lp = 0.0_pReal
|
2014-08-10 15:36:03 +05:30
|
|
|
dLp_dTstar = 0.0_pReal
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_J2_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = 1_pInt
|
|
|
|
call constitutive_j2_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = constitutive_phenopowerlaw_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))
|
|
|
|
call constitutive_phenopowerlaw_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2014-07-02 17:57:39 +05:30
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))
|
|
|
|
call constitutive_nonlocal_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = constitutive_dislotwin_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))
|
|
|
|
call constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2014-08-08 16:34:40 +05:30
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = constitutive_dislokmc_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))
|
|
|
|
call constitutive_dislokmc_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
nSlip = constitutive_titanmod_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))
|
|
|
|
call constitutive_titanmod_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v, &
|
|
|
|
constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
constitutive_getSlipDamage(nSlip,Tstar_v,ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
2014-07-02 17:57:39 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
end select
|
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
end subroutine constitutive_LpAndItsTangent
|
2009-03-06 15:32:36 +05:30
|
|
|
|
|
|
|
|
2014-11-01 00:33:08 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the velocity gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine constitutive_LiAndItsTangent(Li, dLi_dTstar, Tstar_v, Lp, ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
phase_damage, &
|
|
|
|
phase_thermal, &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_LdAndItsTangent
|
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_LTAndItsTangent
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola-Kirchhoff stress
|
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Lp !< plastic velocity gradient
|
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
|
|
|
Li !< intermediate velocity gradient
|
|
|
|
real(pReal), intent(out), dimension(9,9) :: &
|
|
|
|
dLi_dTstar !< derivative of Li with respect to Tstar (2nd-order tensor)
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
Li_temp !< intermediate velocity gradient
|
|
|
|
real(pReal), dimension(9,9) :: &
|
|
|
|
dLi_dTstar_temp !< derivative of Li with respect to Tstar (4th-order tensor)
|
|
|
|
|
|
|
|
Li = 0.0_pReal
|
|
|
|
dLi_dTstar = 0.0_pReal
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
call damage_anisoBrittle_LdAndItsTangent(Li_temp, dLi_dTstar_temp, Tstar_v, ipc, ip, el)
|
|
|
|
Li = Li + Li_temp
|
|
|
|
dLi_dTstar = dLi_dTstar + dLi_dTstar_temp
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
call thermal_adiabatic_LTAndItsTangent(Li_temp, dLi_dTstar_temp, Tstar_v, Lp, ipc, ip, el)
|
|
|
|
Li = Li + Li_temp
|
|
|
|
dLi_dTstar = dLi_dTstar + dLi_dTstar_temp
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end subroutine constitutive_LiAndItsTangent
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the intermediate deformation gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
pure function constitutive_getFi(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use math, only: &
|
|
|
|
math_I3, &
|
|
|
|
math_mul33x33
|
|
|
|
use material, only: &
|
|
|
|
phase_damage, &
|
|
|
|
phase_thermal, &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_getFd
|
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_getFT
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
constitutive_getFi !< intermediate deformation gradient
|
|
|
|
|
|
|
|
constitutive_getFi = math_I3
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_getFi = math_mul33x33(constitutive_getFi,damage_anisoBrittle_getFd (ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
constitutive_getFi = math_mul33x33(constitutive_getFi,thermal_adiabatic_getFT (ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getFi
|
|
|
|
|
|
|
|
|
2014-11-03 16:13:36 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the intermediate deformation gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-11-12 22:10:50 +05:30
|
|
|
subroutine constitutive_putFi(Tstar_v, Lp, dt, ipc, ip, el)
|
2014-11-03 16:13:36 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
phase_damage, &
|
2014-11-12 22:10:50 +05:30
|
|
|
phase_thermal, &
|
2014-11-03 16:13:36 +05:30
|
|
|
material_phase, &
|
2014-11-12 22:10:50 +05:30
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID
|
2014-11-03 16:13:36 +05:30
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_putFd
|
2014-11-12 22:10:50 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_putFT
|
2014-11-03 16:13:36 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola-Kirchhoff stress
|
2014-11-12 22:10:50 +05:30
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Lp !< plastic velocity gradient
|
2014-11-03 16:13:36 +05:30
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
dt
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
call damage_anisoBrittle_putFd (Tstar_v, dt, ipc, ip, el)
|
|
|
|
|
|
|
|
end select
|
2014-11-12 22:10:50 +05:30
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
call thermal_adiabatic_putFT (Tstar_v, Lp, dt, ipc, ip, el)
|
|
|
|
|
|
|
|
end select
|
2014-11-03 16:13:36 +05:30
|
|
|
|
|
|
|
end subroutine constitutive_putFi
|
|
|
|
|
|
|
|
|
2014-11-01 00:33:08 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the intermediate deformation gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
pure function constitutive_getFi0(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use math, only: &
|
|
|
|
math_I3, &
|
|
|
|
math_mul33x33
|
|
|
|
use material, only: &
|
|
|
|
phase_damage, &
|
|
|
|
phase_thermal, &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_getFd0
|
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_getFT0
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
constitutive_getFi0 !< intermediate deformation gradient
|
|
|
|
|
|
|
|
constitutive_getFi0 = math_I3
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_getFi0 = math_mul33x33(constitutive_getFi0,damage_anisoBrittle_getFd0 (ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
constitutive_getFi0 = math_mul33x33(constitutive_getFi0,thermal_adiabatic_getFT0 (ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getFi0
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the intermediate deformation gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
pure function constitutive_getPartionedFi0(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use math, only: &
|
|
|
|
math_I3, &
|
|
|
|
math_mul33x33
|
|
|
|
use material, only: &
|
|
|
|
phase_damage, &
|
|
|
|
phase_thermal, &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_getPartionedFd0
|
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_getPartionedFT0
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
constitutive_getPartionedFi0 !< intermediate deformation gradient
|
|
|
|
|
|
|
|
constitutive_getPartionedFi0 = math_I3
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_getPartionedFi0 = math_mul33x33(constitutive_getPartionedFi0, &
|
|
|
|
damage_anisoBrittle_getPartionedFd0(ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
constitutive_getPartionedFi0 = math_mul33x33(constitutive_getPartionedFi0, &
|
|
|
|
thermal_adiabatic_getPartionedFT0(ipc, ip, el))
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getPartionedFi0
|
|
|
|
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to
|
2013-10-14 16:24:45 +05:30
|
|
|
!> the elastic deformation gradient depending on the selected elastic law (so far no case switch
|
|
|
|
!! because only hooke is implemented
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-07-02 17:57:39 +05:30
|
|
|
subroutine constitutive_TandItsTangent(T, dT_dFe, Fe, ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Fe !< elastic deformation gradient
|
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
|
|
|
T !< 2nd Piola-Kirchhoff stress tensor
|
|
|
|
real(pReal), intent(out), dimension(3,3,3,3) :: &
|
|
|
|
dT_dFe !< derivative of 2nd P-K stress with respect to elastic deformation gradient
|
|
|
|
|
2013-10-14 16:24:45 +05:30
|
|
|
call constitutive_hooke_TandItsTangent(T, dT_dFe, Fe, ipc, ip, el)
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
end subroutine constitutive_TandItsTangent
|
2012-03-15 15:21:33 +05:30
|
|
|
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to
|
2013-10-14 20:05:41 +05:30
|
|
|
!> the elastic deformation gradient using hookes law
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-07-02 17:57:39 +05:30
|
|
|
subroutine constitutive_hooke_TandItsTangent(T, dT_dFe, Fe, ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use math, only : &
|
|
|
|
math_mul3x3, &
|
|
|
|
math_mul33x33, &
|
|
|
|
math_mul3333xx33, &
|
|
|
|
math_Mandel66to3333, &
|
|
|
|
math_transpose33, &
|
2014-08-09 03:59:38 +05:30
|
|
|
math_trace33, &
|
|
|
|
math_I3
|
2014-07-02 17:57:39 +05:30
|
|
|
use material, only: &
|
2014-08-08 19:32:20 +05:30
|
|
|
mappingConstitutive
|
2014-06-25 04:52:52 +05:30
|
|
|
use lattice, only: &
|
2014-07-02 17:57:39 +05:30
|
|
|
lattice_referenceTemperature, &
|
|
|
|
lattice_thermalExpansion33
|
2012-11-07 21:13:29 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Fe !< elastic deformation gradient
|
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
|
|
|
T !< 2nd Piola-Kirchhoff stress tensor
|
|
|
|
real(pReal), intent(out), dimension(3,3,3,3) :: &
|
|
|
|
dT_dFe !< dT/dFe
|
|
|
|
|
2014-01-22 14:08:13 +05:30
|
|
|
integer(pInt) :: i, j, k, l
|
2014-08-09 03:59:38 +05:30
|
|
|
real(pReal), dimension(3,3,3,3) :: C
|
|
|
|
|
2014-10-27 21:03:35 +05:30
|
|
|
C = math_Mandel66to3333(constitutive_damagedC(ipc,ip,el))
|
2014-10-20 21:13:28 +05:30
|
|
|
T = math_mul3333xx33(C,0.5_pReal*(math_mul33x33(math_transpose33(Fe),Fe)-math_I3))
|
2014-08-10 16:15:07 +05:30
|
|
|
|
2014-09-23 02:08:19 +05:30
|
|
|
dT_dFe = 0.0_pReal
|
|
|
|
forall (i=1_pInt:3_pInt, j=1_pInt:3_pInt, k=1_pInt:3_pInt, l=1_pInt:3_pInt) &
|
|
|
|
dT_dFe(i,j,k,l) = sum(C(i,j,l,1:3)*Fe(k,1:3)) ! dT*_ij/dFe_kl
|
2014-08-10 16:15:07 +05:30
|
|
|
|
2012-03-15 15:21:33 +05:30
|
|
|
end subroutine constitutive_hooke_TandItsTangent
|
|
|
|
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the rate of change of microstructure
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-10 17:58:57 +05:30
|
|
|
subroutine constitutive_collectDotState(Tstar_v, Lp, FeArray, FpArray, subdt, subfracArray,&
|
2014-03-13 12:13:49 +05:30
|
|
|
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_cumDotStateCalls, &
|
|
|
|
debug_cumDotStateTicks, &
|
|
|
|
debug_level, &
|
|
|
|
debug_constitutive, &
|
|
|
|
debug_levelBasic
|
|
|
|
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: &
|
|
|
|
phase_plasticity, &
|
2014-09-23 16:08:20 +05:30
|
|
|
phase_damage, &
|
2014-10-11 02:25:09 +05:30
|
|
|
phase_vacancy, &
|
2013-01-16 15:44:57 +05:30
|
|
|
material_phase, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
2014-10-09 19:38:32 +05:30
|
|
|
PLASTICITY_none_ID, &
|
|
|
|
PLASTICITY_j2_ID, &
|
|
|
|
PLASTICITY_phenopowerlaw_ID, &
|
|
|
|
PLASTICITY_dislotwin_ID, &
|
|
|
|
PLASTICITY_dislokmc_ID, &
|
|
|
|
PLASTICITY_titanmod_ID, &
|
|
|
|
PLASTICITY_nonlocal_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
2014-10-11 02:25:09 +05:30
|
|
|
LOCAL_VACANCY_generation_ID
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_j2, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_j2_dotState
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_phenopowerlaw, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_phenopowerlaw_dotState
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_dislotwin, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_dislotwin_dotState
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_dotState
|
2013-11-27 13:34:05 +05:30
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_dotState
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_nonlocal, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_nonlocal_dotState
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_dotState
|
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_dotState
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_dotState
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile, only: &
|
|
|
|
damage_anisoDuctile_dotState
|
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_dotState
|
2014-10-11 02:25:09 +05:30
|
|
|
use vacancy_generation, only: &
|
|
|
|
vacancy_generation_dotState
|
2014-06-14 02:23:17 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
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
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola Kirchhoff stress tensor (Mandel)
|
2014-09-27 02:19:25 +05:30
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Lp !< plastic velocity gradient
|
2013-01-16 15:44:57 +05:30
|
|
|
integer(pLongInt) :: &
|
|
|
|
tick, tock, &
|
|
|
|
tickrate, &
|
|
|
|
maxticks
|
2014-10-11 15:15:30 +05:30
|
|
|
real(pReal), dimension(:), allocatable :: &
|
|
|
|
accumulatedSlip
|
|
|
|
integer(pInt) :: &
|
|
|
|
nSlip
|
2012-07-03 16:46:38 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
if (iand(debug_level(debug_constitutive), debug_levelBasic) /= 0_pInt) &
|
|
|
|
call system_clock(count=tick,count_rate=tickrate,count_max=maxticks)
|
2009-12-15 13:50:31 +05:30
|
|
|
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_J2_ID)
|
2014-08-10 16:17:12 +05:30
|
|
|
call constitutive_j2_dotState (Tstar_v,ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
2014-08-10 16:17:12 +05:30
|
|
|
call constitutive_phenopowerlaw_dotState(Tstar_v,ipc,ip,el)
|
2014-06-03 19:16:42 +05:30
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_dislotwin_dotState (Tstar_v,constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2014-08-08 16:34:40 +05:30
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_dislokmc_dotState (Tstar_v,constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_titanmod_dotState (Tstar_v,constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
2014-10-10 17:58:57 +05:30
|
|
|
call constitutive_nonlocal_dotState (Tstar_v,FeArray,FpArray,constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
subdt,subfracArray,ip,el)
|
2013-01-16 15:44:57 +05:30
|
|
|
end select
|
|
|
|
|
2014-09-23 16:08:20 +05:30
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
2014-11-12 04:16:41 +05:30
|
|
|
call damage_isoBrittle_dotState(constitutive_homogenizedC(ipc,ip,el), &
|
|
|
|
FeArray(1:3,1:3,ipc,ip,el), ipc, ip, el)
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
2014-11-12 04:16:41 +05:30
|
|
|
call constitutive_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
call damage_isoDuctile_dotState(nSlip,accumulatedSlip,ipc, ip, el)
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
call damage_anisoBrittle_dotState(Tstar_v, ipc, ip, el)
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
2014-11-06 23:22:43 +05:30
|
|
|
call constitutive_getAccumulatedSlip(nSlip,accumulatedSlip,ipc,ip,el)
|
|
|
|
call damage_anisoDuctile_dotState(nSlip, accumulatedSlip, ipc, ip, el)
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
call damage_gurson_dotState(Tstar_v, Lp, ipc, ip, el)
|
2014-09-23 16:08:20 +05:30
|
|
|
end select
|
|
|
|
|
2014-10-11 02:25:09 +05:30
|
|
|
select case (phase_vacancy(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
2014-11-05 23:11:08 +05:30
|
|
|
call constitutive_getAccumulatedSlip(nSlip,accumulatedSlip,ipc,ip,el)
|
2014-10-11 15:15:30 +05:30
|
|
|
call vacancy_generation_dotState(nSlip,accumulatedSlip,Tstar_v,constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
ipc, ip, el)
|
2014-10-11 02:25:09 +05:30
|
|
|
end select
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
if (iand(debug_level(debug_constitutive), debug_levelBasic) /= 0_pInt) then
|
|
|
|
call system_clock(count=tock,count_rate=tickrate,count_max=maxticks)
|
|
|
|
!$OMP CRITICAL (debugTimingDotState)
|
|
|
|
debug_cumDotStateCalls = debug_cumDotStateCalls + 1_pInt
|
|
|
|
debug_cumDotStateTicks = debug_cumDotStateTicks + tock-tick
|
|
|
|
!$OMP FLUSH (debug_cumDotStateTicks)
|
|
|
|
if (tock < tick) debug_cumDotStateTicks = debug_cumDotStateTicks + maxticks
|
|
|
|
!$OMP END CRITICAL (debugTimingDotState)
|
|
|
|
endif
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-06-17 12:24:49 +05:30
|
|
|
!> @brief for constitutive models having an instantaneous change of state (so far, only nonlocal)
|
|
|
|
!> will return false if delta state is not needed/supported by the constitutive model
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-06-17 12:24:49 +05:30
|
|
|
logical function constitutive_collectDeltaState(Tstar_v, 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_cumDeltaStateCalls, &
|
|
|
|
debug_cumDeltaStateTicks, &
|
|
|
|
debug_level, &
|
|
|
|
debug_constitutive, &
|
|
|
|
debug_levelBasic
|
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
2013-11-27 13:34:05 +05:30
|
|
|
material_phase, &
|
2014-06-14 02:23:17 +05:30
|
|
|
plasticState, &
|
2014-06-23 00:28:29 +05:30
|
|
|
mappingConstitutive, &
|
2014-06-14 02:23:17 +05:30
|
|
|
PLASTICITY_NONLOCAL_ID
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_nonlocal, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_nonlocal_deltaState
|
2013-01-16 15:44:57 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2013-01-16 15:44:57 +05:30
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola-Kirchhoff stress
|
|
|
|
integer(pLongInt) :: &
|
|
|
|
tick, tock, &
|
|
|
|
tickrate, &
|
|
|
|
maxticks
|
|
|
|
|
|
|
|
if (iand(debug_level(debug_constitutive), debug_levelBasic) /= 0_pInt) &
|
|
|
|
call system_clock(count=tick,count_rate=tickrate,count_max=maxticks)
|
|
|
|
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2012-05-16 20:13:26 +05:30
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
2014-06-17 12:24:49 +05:30
|
|
|
constitutive_collectDeltaState = .true.
|
2014-06-26 19:23:12 +05:30
|
|
|
call constitutive_nonlocal_deltaState(Tstar_v,ip,el)
|
2013-10-14 16:24:45 +05:30
|
|
|
case default
|
2014-06-17 12:24:49 +05:30
|
|
|
constitutive_collectDeltaState = .false.
|
2014-07-02 17:57:39 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
end select
|
2012-05-16 20:13:26 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
if (iand(debug_level(debug_constitutive), debug_levelBasic) /= 0_pInt) then
|
|
|
|
call system_clock(count=tock,count_rate=tickrate,count_max=maxticks)
|
|
|
|
!$OMP CRITICAL (debugTimingDeltaState)
|
|
|
|
debug_cumDeltaStateCalls = debug_cumDeltaStateCalls + 1_pInt
|
|
|
|
debug_cumDeltaStateTicks = debug_cumDeltaStateTicks + tock-tick
|
|
|
|
!$OMP FLUSH (debug_cumDeltaStateTicks)
|
|
|
|
if (tock < tick) debug_cumDeltaStateTicks = debug_cumDeltaStateTicks + maxticks
|
|
|
|
!$OMP END CRITICAL (debugTimingDeltaState)
|
|
|
|
endif
|
2012-05-16 20:13:26 +05:30
|
|
|
|
2014-06-17 12:24:49 +05:30
|
|
|
end function constitutive_collectDeltaState
|
2014-09-10 23:56:12 +05:30
|
|
|
|
2014-09-22 23:45:19 +05:30
|
|
|
|
2014-09-05 22:01:27 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-15 20:32:29 +05:30
|
|
|
!> @brief Returns the local(regularised) damage
|
2014-09-05 22:01:27 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_getLocalDamage(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_none_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
2014-09-05 22:01:27 +05:30
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_getLocalDamage
|
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_getLocalDamage
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_getLocalDamage
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile, only: &
|
|
|
|
damage_anisoDuctile_getLocalDamage
|
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_getLocalDamage
|
2014-09-05 22:01:27 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal) :: constitutive_getLocalDamage
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_none_ID)
|
|
|
|
constitutive_getLocalDamage = 1.0_pReal
|
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
constitutive_getLocalDamage = damage_isoBrittle_getLocalDamage(ipc, ip, el)
|
2014-09-26 16:08:13 +05:30
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
constitutive_getLocalDamage = damage_isoDuctile_getLocalDamage(ipc, ip, el)
|
2014-10-10 18:12:12 +05:30
|
|
|
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_getLocalDamage = damage_anisoBrittle_getLocalDamage(ipc, ip, el)
|
|
|
|
|
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
|
|
|
|
|
|
|
constitutive_getLocalDamage = damage_anisoDuctile_getLocalDamage(ipc, ip, el)
|
|
|
|
|
2014-10-10 18:12:12 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
2014-10-28 06:48:10 +05:30
|
|
|
constitutive_getLocalDamage = damage_gurson_getLocalDamage(ipc, ip, el)
|
2014-10-15 20:32:29 +05:30
|
|
|
|
2014-09-05 22:01:27 +05:30
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getLocalDamage
|
|
|
|
|
2014-09-22 23:45:19 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Returns the local(unregularised) damage
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine constitutive_putLocalDamage(ipc, ip, el, localDamage)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
2014-09-22 23:45:19 +05:30
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_putLocalDamage
|
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_putLocalDamage
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_putLocalDamage
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile, only: &
|
|
|
|
damage_anisoDuctile_putLocalDamage
|
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_putLocalDamage
|
2014-09-22 23:45:19 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
localDamage
|
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
call damage_isoBrittle_putLocalDamage(ipc, ip, el, localDamage)
|
2014-10-06 22:31:39 +05:30
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
call damage_isoDuctile_putLocalDamage(ipc, ip, el, localDamage)
|
2014-10-10 18:12:12 +05:30
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
call damage_anisoBrittle_putLocalDamage(ipc, ip, el, localDamage)
|
2014-09-22 23:45:19 +05:30
|
|
|
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
|
|
|
call damage_anisoDuctile_putLocalDamage(ipc, ip, el, localDamage)
|
|
|
|
|
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
call damage_gurson_putLocalDamage(ipc, ip, el, localDamage)
|
|
|
|
|
2014-09-22 23:45:19 +05:30
|
|
|
end select
|
|
|
|
|
|
|
|
end subroutine constitutive_putLocalDamage
|
|
|
|
|
2014-09-10 23:56:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns nonlocal (regularised) damage
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-09-26 23:37:48 +05:30
|
|
|
function constitutive_getDamage(ipc, ip, el)
|
2014-09-05 22:01:27 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
2014-10-28 06:48:10 +05:30
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_none_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-28 06:48:10 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_getDamage
|
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_getDamage
|
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_getDamage
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile, only: &
|
|
|
|
damage_anisoDuctile_getDamage
|
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_getDamage
|
2014-09-05 22:01:27 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2014-09-26 23:37:48 +05:30
|
|
|
real(pReal) :: constitutive_getDamage
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-10-28 06:48:10 +05:30
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_DAMAGE_none_ID)
|
|
|
|
constitutive_getDamage = 1.0_pReal
|
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
constitutive_getDamage = damage_isoBrittle_getDamage(ipc, ip, el)
|
2014-10-28 06:48:10 +05:30
|
|
|
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
constitutive_getDamage = damage_isoDuctile_getDamage(ipc, ip, el)
|
2014-10-28 06:48:10 +05:30
|
|
|
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_getDamage = damage_anisoBrittle_getDamage(ipc, ip, el)
|
|
|
|
|
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
|
|
|
constitutive_getDamage = damage_anisoDuctile_getDamage(ipc, ip, el)
|
|
|
|
|
2014-10-28 06:48:10 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
constitutive_getDamage = damage_gurson_getDamage(ipc, ip, el)
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-09-23 17:52:34 +05:30
|
|
|
end select
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-09-26 23:37:48 +05:30
|
|
|
end function constitutive_getDamage
|
2014-10-11 15:39:36 +05:30
|
|
|
|
2014-10-28 08:12:25 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Returns the damage on each slip system
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_getSlipDamage(nSlip, Tstar_v, ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
2014-11-05 23:11:08 +05:30
|
|
|
LOCAL_DAMAGE_anisoDuctile_ID, &
|
2014-10-28 08:12:25 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_getSlipDamage
|
2014-11-05 23:11:08 +05:30
|
|
|
use damage_anisoDuctile, only: &
|
|
|
|
damage_anisoDuctile_getSlipDamage
|
2014-10-28 08:12:25 +05:30
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_getSlipDamage
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
nSlip, &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in), dimension(6) :: &
|
|
|
|
Tstar_v !< 2nd Piola-Kirchhoff stress
|
2014-11-05 23:11:08 +05:30
|
|
|
real(pReal) :: &
|
|
|
|
constitutive_getSlipDamage(nSlip)
|
2014-10-28 08:12:25 +05:30
|
|
|
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
constitutive_getSlipDamage = damage_isoDuctile_getSlipDamage(ipc, ip, el)
|
2014-10-28 08:12:25 +05:30
|
|
|
|
2014-11-05 23:11:08 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoDuctile_ID)
|
2014-11-06 23:22:43 +05:30
|
|
|
constitutive_getSlipDamage = damage_anisoDuctile_getSlipDamage(ipc, ip, el)
|
2014-11-05 23:11:08 +05:30
|
|
|
|
2014-10-28 08:12:25 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
constitutive_getSlipDamage = damage_gurson_getSlipDamage(Tstar_v, ipc, ip, el)
|
|
|
|
|
|
|
|
case default
|
|
|
|
constitutive_getSlipDamage = 1.0_pReal
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getSlipDamage
|
|
|
|
|
2014-10-11 15:39:36 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns damage diffusion tensor
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_getDamageDiffusion33(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use lattice, only: &
|
|
|
|
lattice_DamageDiffusion33
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_DAMAGE_none_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-10-11 15:39:36 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
|
|
|
phase_damage
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_getDamageDiffusion33
|
2014-10-11 15:39:36 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
constitutive_getDamageDiffusion33
|
|
|
|
|
|
|
|
constitutive_getDamageDiffusion33 = lattice_DamageDiffusion33(1:3,1:3,material_phase(ipc,ip,el))
|
|
|
|
select case(phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
constitutive_getDamageDiffusion33 = damage_isoBrittle_getDamageDiffusion33(ipc, ip, el)
|
2014-10-11 15:39:36 +05:30
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getDamageDiffusion33
|
|
|
|
|
2014-09-10 23:56:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns local (unregularised) temperature
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-09-26 23:37:48 +05:30
|
|
|
function constitutive_getAdiabaticTemperature(ipc, ip, el)
|
2014-09-05 22:01:27 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
2014-10-09 19:38:32 +05:30
|
|
|
LOCAL_THERMAL_isothermal_ID, &
|
|
|
|
LOCAL_THERMAL_adiabatic_ID, &
|
2014-10-13 23:24:27 +05:30
|
|
|
phase_thermal, &
|
|
|
|
phase_thermalInstance
|
|
|
|
use thermal_isothermal, only: &
|
|
|
|
thermal_isothermal_temperature
|
2014-10-09 19:38:32 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_getTemperature
|
2014-09-05 22:01:27 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2014-09-26 23:37:48 +05:30
|
|
|
real(pReal) :: constitutive_getAdiabaticTemperature
|
2014-09-05 22:01:27 +05:30
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
2014-10-10 14:10:59 +05:30
|
|
|
case (LOCAL_THERMAL_isothermal_ID)
|
2014-10-13 23:24:27 +05:30
|
|
|
constitutive_getAdiabaticTemperature = &
|
|
|
|
thermal_isothermal_temperature(phase_thermalInstance(material_phase(ipc,ip,el)))
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-10-09 19:38:32 +05:30
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
constitutive_getAdiabaticTemperature = thermal_adiabatic_getTemperature(ipc, ip, el)
|
2014-09-05 22:01:27 +05:30
|
|
|
end select
|
|
|
|
|
2014-09-26 23:37:48 +05:30
|
|
|
end function constitutive_getAdiabaticTemperature
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-09-22 23:45:19 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-28 16:19:12 +05:30
|
|
|
!> @brief assigns the local/nonlocal value of temperature to local thermal state
|
2014-09-22 23:45:19 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-09-26 23:37:48 +05:30
|
|
|
subroutine constitutive_putAdiabaticTemperature(ipc, ip, el, localTemperature)
|
2014-09-22 23:45:19 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
2014-10-09 19:38:32 +05:30
|
|
|
LOCAL_THERMAL_adiabatic_ID, &
|
2014-09-22 23:45:19 +05:30
|
|
|
phase_thermal
|
2014-10-09 19:38:32 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_putTemperature
|
2014-09-22 23:45:19 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
localTemperature
|
|
|
|
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
2014-10-09 19:38:32 +05:30
|
|
|
case (LOCAL_THERMAL_adiabatic_ID)
|
|
|
|
call thermal_adiabatic_putTemperature(ipc, ip, el, localTemperature)
|
2014-09-22 23:45:19 +05:30
|
|
|
|
|
|
|
end select
|
|
|
|
|
2014-09-26 23:37:48 +05:30
|
|
|
end subroutine constitutive_putAdiabaticTemperature
|
2014-09-22 23:45:19 +05:30
|
|
|
|
2014-09-10 23:56:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns nonlocal (regularised) temperature
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-09-26 23:37:48 +05:30
|
|
|
function constitutive_getTemperature(ipc, ip, el)
|
2014-09-05 22:01:27 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
mappingHomogenization, &
|
2014-09-10 20:42:14 +05:30
|
|
|
material_phase, &
|
2014-09-05 22:01:27 +05:30
|
|
|
fieldThermal, &
|
|
|
|
field_thermal_type, &
|
2014-10-09 19:38:32 +05:30
|
|
|
FIELD_THERMAL_local_ID, &
|
|
|
|
FIELD_THERMAL_nonlocal_ID, &
|
2014-09-05 22:01:27 +05:30
|
|
|
material_homog
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2014-09-26 23:37:48 +05:30
|
|
|
real(pReal) :: constitutive_getTemperature
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-10-28 06:48:10 +05:30
|
|
|
select case(field_thermal_type(material_homog(ip,el)))
|
|
|
|
case (FIELD_THERMAL_local_ID)
|
|
|
|
constitutive_getTemperature = constitutive_getAdiabaticTemperature(ipc, ip, el)
|
|
|
|
|
|
|
|
case (FIELD_THERMAL_nonlocal_ID)
|
|
|
|
constitutive_getTemperature = fieldThermal(material_homog(ip,el))% &
|
|
|
|
field(1,mappingHomogenization(1,ip,el)) ! Taylor type
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-10-28 06:48:10 +05:30
|
|
|
end select
|
2014-09-05 22:01:27 +05:30
|
|
|
|
2014-09-26 23:37:48 +05:30
|
|
|
end function constitutive_getTemperature
|
2014-10-10 22:04:51 +05:30
|
|
|
|
2014-10-11 02:25:09 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns local vacancy concentration
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_getLocalVacancyConcentration(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_VACANCY_constant_ID, &
|
|
|
|
LOCAL_VACANCY_generation_ID, &
|
|
|
|
phase_vacancy
|
|
|
|
use vacancy_generation, only: &
|
|
|
|
vacancy_generation_getConcentration
|
|
|
|
use lattice, only: &
|
|
|
|
lattice_equilibriumVacancyConcentration
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal) :: constitutive_getLocalVacancyConcentration
|
|
|
|
|
|
|
|
select case (phase_vacancy(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_VACANCY_constant_ID)
|
|
|
|
constitutive_getLocalVacancyConcentration = &
|
|
|
|
lattice_equilibriumVacancyConcentration(material_phase(ipc,ip,el))
|
|
|
|
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
|
|
|
constitutive_getLocalVacancyConcentration = vacancy_generation_getConcentration(ipc, ip, el)
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getLocalVacancyConcentration
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Puts local vacancy concentration
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine constitutive_putLocalVacancyConcentration(ipc, ip, el, localVacancyConcentration)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_VACANCY_generation_ID, &
|
|
|
|
phase_vacancy
|
|
|
|
use vacancy_generation, only: &
|
|
|
|
vacancy_generation_putConcentration
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
localVacancyConcentration
|
|
|
|
|
|
|
|
select case (phase_vacancy(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
|
|
|
call vacancy_generation_putConcentration(ipc, ip, el, localVacancyConcentration)
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end subroutine constitutive_putLocalVacancyConcentration
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns nonlocal vacancy concentration
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function constitutive_getVacancyConcentration(ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use material, only: &
|
|
|
|
mappingHomogenization, &
|
|
|
|
material_phase, &
|
|
|
|
fieldVacancy, &
|
|
|
|
field_vacancy_type, &
|
|
|
|
FIELD_VACANCY_local_ID, &
|
|
|
|
FIELD_VACANCY_nonlocal_ID, &
|
|
|
|
material_homog
|
|
|
|
use lattice, only: &
|
|
|
|
lattice_equilibriumVacancyConcentration
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal) :: constitutive_getVacancyConcentration
|
|
|
|
|
|
|
|
select case(field_vacancy_type(material_homog(ip,el)))
|
|
|
|
case (FIELD_VACANCY_local_ID)
|
|
|
|
constitutive_getVacancyConcentration = constitutive_getLocalVacancyConcentration(ipc, ip, el)
|
|
|
|
|
|
|
|
case (FIELD_VACANCY_nonlocal_ID)
|
|
|
|
constitutive_getVacancyConcentration = fieldVacancy(material_homog(ip,el))% &
|
|
|
|
field(1,mappingHomogenization(1,ip,el)) ! Taylor type
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getVacancyConcentration
|
|
|
|
|
2014-10-11 16:09:44 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns vacancy diffusion tensor
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-11-05 23:11:08 +05:30
|
|
|
function constitutive_getVacancyDiffusion33(ipc, ip, el)
|
2014-10-11 16:09:44 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
|
|
|
use lattice, only: &
|
|
|
|
lattice_VacancyDiffusion33
|
|
|
|
use material, only: &
|
|
|
|
material_phase, &
|
|
|
|
LOCAL_VACANCY_generation_ID, &
|
|
|
|
phase_vacancy
|
|
|
|
use vacancy_generation, only: &
|
|
|
|
vacancy_generation_getVacancyDiffusion33
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
constitutive_getVacancyDiffusion33
|
|
|
|
real(pReal), dimension(:), allocatable :: &
|
|
|
|
accumulatedSlip
|
|
|
|
integer(pInt) :: &
|
|
|
|
nSlip
|
|
|
|
|
|
|
|
constitutive_getVacancyDiffusion33 = lattice_VacancyDiffusion33(1:3,1:3,material_phase(ipc,ip,el))
|
|
|
|
select case(phase_vacancy(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
2014-11-05 23:11:08 +05:30
|
|
|
call constitutive_getAccumulatedSlip(nSlip,accumulatedSlip,ipc,ip,el)
|
2014-10-11 16:09:44 +05:30
|
|
|
constitutive_getVacancyDiffusion33 = &
|
|
|
|
vacancy_generation_getVacancyDiffusion33(nSlip,accumulatedSlip,constitutive_getTemperature(ipc,ip,el), &
|
|
|
|
ipc,ip,el)
|
|
|
|
|
|
|
|
end select
|
|
|
|
|
|
|
|
end function constitutive_getVacancyDiffusion33
|
|
|
|
|
2014-10-10 22:04:51 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns accumulated slip on each system defined
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-11-05 23:11:08 +05:30
|
|
|
subroutine constitutive_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
2014-10-10 22:04:51 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal, &
|
|
|
|
pInt
|
|
|
|
use math, only: &
|
|
|
|
math_mul33xx33, &
|
2014-10-11 15:15:30 +05:30
|
|
|
math_equivStrain33, &
|
2014-10-10 22:04:51 +05:30
|
|
|
math_I3
|
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
|
|
|
material_phase, &
|
|
|
|
PLASTICITY_none_ID, &
|
|
|
|
PLASTICITY_j2_ID, &
|
|
|
|
PLASTICITY_phenopowerlaw_ID, &
|
|
|
|
PLASTICITY_dislotwin_ID, &
|
|
|
|
PLASTICITY_dislokmc_ID, &
|
|
|
|
PLASTICITY_titanmod_ID, &
|
|
|
|
PLASTICITY_nonlocal_ID
|
2014-11-05 23:11:08 +05:30
|
|
|
use constitutive_J2, only: &
|
|
|
|
constitutive_J2_getAccumulatedSlip
|
2014-10-10 22:04:51 +05:30
|
|
|
use constitutive_phenopowerlaw, only: &
|
|
|
|
constitutive_phenopowerlaw_getAccumulatedSlip
|
|
|
|
use constitutive_dislotwin, only: &
|
|
|
|
constitutive_dislotwin_getAccumulatedSlip
|
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_getAccumulatedSlip
|
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_getAccumulatedSlip
|
|
|
|
use constitutive_nonlocal, only: &
|
|
|
|
constitutive_nonlocal_getAccumulatedSlip
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
real(pReal), dimension(:), allocatable :: &
|
|
|
|
accumulatedSlip
|
|
|
|
integer(pInt) :: &
|
|
|
|
nSlip
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
|
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
|
|
|
case (PLASTICITY_none_ID)
|
|
|
|
nSlip = 0_pInt
|
|
|
|
allocate(accumulatedSlip(nSlip))
|
|
|
|
case (PLASTICITY_J2_ID)
|
2014-11-05 23:11:08 +05:30
|
|
|
call constitutive_J2_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
2014-10-10 22:04:51 +05:30
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
|
|
|
call constitutive_phenopowerlaw_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
|
|
|
call constitutive_dislotwin_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
|
|
|
call constitutive_dislokmc_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
|
|
|
call constitutive_titanmod_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
|
|
|
call constitutive_nonlocal_getAccumulatedSlip(nSlip,accumulatedSlip,ipc, ip, el)
|
|
|
|
end select
|
|
|
|
|
|
|
|
end subroutine constitutive_getAccumulatedSlip
|
|
|
|
|
2014-10-11 15:15:30 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns accumulated slip rates on each system defined
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine constitutive_getSlipRate(nSlip,slipRate,Lp,ipc, ip, el)
|
|
|
|
use prec, only: &
|
|
|
|
pReal, &
|
|
|
|
pInt
|
|
|
|
use math, only: &
|
|
|
|
math_mul33xx33, &
|
|
|
|
math_equivStrain33, &
|
|
|
|
math_I3
|
|
|
|
use material, only: &
|
|
|
|
phase_plasticity, &
|
|
|
|
material_phase, &
|
|
|
|
PLASTICITY_none_ID, &
|
|
|
|
PLASTICITY_j2_ID, &
|
|
|
|
PLASTICITY_phenopowerlaw_ID, &
|
|
|
|
PLASTICITY_dislotwin_ID, &
|
|
|
|
PLASTICITY_dislokmc_ID, &
|
|
|
|
PLASTICITY_titanmod_ID, &
|
|
|
|
PLASTICITY_nonlocal_ID
|
|
|
|
use constitutive_phenopowerlaw, only: &
|
|
|
|
constitutive_phenopowerlaw_getSlipRate
|
|
|
|
use constitutive_dislotwin, only: &
|
|
|
|
constitutive_dislotwin_getSlipRate
|
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_getSlipRate
|
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_getSlipRate
|
|
|
|
use constitutive_nonlocal, only: &
|
|
|
|
constitutive_nonlocal_getSlipRate
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
real(pReal), dimension(:), allocatable :: &
|
|
|
|
slipRate
|
|
|
|
integer(pInt) :: &
|
|
|
|
nSlip
|
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
Lp !< plastic velocity gradient
|
|
|
|
integer(pInt), intent(in) :: &
|
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
|
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
|
|
|
case (PLASTICITY_none_ID)
|
|
|
|
nSlip = 0_pInt
|
|
|
|
allocate(slipRate(nSlip))
|
|
|
|
case (PLASTICITY_J2_ID)
|
|
|
|
nSlip = 1_pInt
|
|
|
|
allocate(slipRate(nSlip))
|
|
|
|
slipRate(1) = math_equivStrain33(Lp)
|
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
|
|
|
call constitutive_phenopowerlaw_getSlipRate(nSlip,slipRate,ipc, ip, el)
|
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
|
|
|
call constitutive_dislotwin_getSlipRate(nSlip,slipRate,ipc, ip, el)
|
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
|
|
|
call constitutive_dislokmc_getSlipRate(nSlip,slipRate,ipc, ip, el)
|
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
|
|
|
call constitutive_titanmod_getSlipRate(nSlip,slipRate,ipc, ip, el)
|
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
|
|
|
call constitutive_nonlocal_getSlipRate(nSlip,slipRate,ipc, ip, el)
|
|
|
|
end select
|
|
|
|
|
|
|
|
end subroutine constitutive_getSlipRate
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns array of constitutive results
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2014-10-10 17:58:57 +05:30
|
|
|
function constitutive_postResults(Tstar_v, FeArray, ipc, ip, el)
|
2014-07-02 17:57:39 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-01-16 15:44:57 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_NcpElems, &
|
|
|
|
mesh_maxNips
|
|
|
|
use material, only: &
|
2014-05-22 20:54:12 +05:30
|
|
|
plasticState, &
|
2014-09-23 16:08:20 +05:30
|
|
|
damageState, &
|
|
|
|
thermalState, &
|
2014-10-11 02:25:09 +05:30
|
|
|
vacancyState, &
|
2013-01-16 15:44:57 +05:30
|
|
|
phase_plasticity, &
|
2014-09-23 16:08:20 +05:30
|
|
|
phase_damage, &
|
|
|
|
phase_thermal, &
|
2014-10-11 02:25:09 +05:30
|
|
|
phase_vacancy, &
|
2013-01-16 15:44:57 +05:30
|
|
|
material_phase, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
|
|
|
PLASTICITY_NONE_ID, &
|
|
|
|
PLASTICITY_J2_ID, &
|
|
|
|
PLASTICITY_PHENOPOWERLAW_ID, &
|
|
|
|
PLASTICITY_DISLOTWIN_ID, &
|
2014-08-08 16:34:40 +05:30
|
|
|
PLASTICITY_DISLOKMC_ID, &
|
2013-11-27 13:34:05 +05:30
|
|
|
PLASTICITY_TITANMOD_ID, &
|
2014-09-23 16:08:20 +05:30
|
|
|
PLASTICITY_NONLOCAL_ID, &
|
2014-10-28 16:19:12 +05:30
|
|
|
LOCAL_DAMAGE_isoBrittle_ID, &
|
|
|
|
LOCAL_DAMAGE_isoDuctile_ID, &
|
|
|
|
LOCAL_DAMAGE_anisoBrittle_ID, &
|
2014-10-10 18:12:12 +05:30
|
|
|
LOCAL_DAMAGE_gurson_ID, &
|
2014-10-11 02:25:09 +05:30
|
|
|
LOCAL_THERMAL_ADIABATIC_ID, &
|
|
|
|
LOCAL_VACANCY_generation_ID
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_j2, only: &
|
2014-04-15 15:39:20 +05:30
|
|
|
#ifdef HDF
|
|
|
|
constitutive_j2_postResults2,&
|
|
|
|
#endif
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_j2_postResults
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_phenopowerlaw, only: &
|
2014-05-22 20:54:12 +05:30
|
|
|
constitutive_phenopowerlaw_postResults
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_dislotwin, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_dislotwin_postResults
|
2014-08-08 16:34:40 +05:30
|
|
|
use constitutive_dislokmc, only: &
|
|
|
|
constitutive_dislokmc_postResults
|
2013-11-27 13:34:05 +05:30
|
|
|
use constitutive_titanmod, only: &
|
|
|
|
constitutive_titanmod_postResults
|
2013-01-16 15:44:57 +05:30
|
|
|
use constitutive_nonlocal, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
constitutive_nonlocal_postResults
|
2014-09-23 16:08:20 +05:30
|
|
|
#ifdef multiphysicsOut
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_isoBrittle, only: &
|
|
|
|
damage_isoBrittle_postResults
|
|
|
|
use damage_isoDuctile, only: &
|
|
|
|
damage_isoDuctile_postResults
|
2014-10-10 18:12:12 +05:30
|
|
|
use damage_gurson, only: &
|
|
|
|
damage_gurson_postResults
|
2014-10-28 16:19:12 +05:30
|
|
|
use damage_anisoBrittle, only: &
|
|
|
|
damage_anisoBrittle_postResults
|
2014-10-09 19:38:32 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_postResults
|
2014-10-11 02:25:09 +05:30
|
|
|
use vacancy_generation, only: &
|
|
|
|
vacancy_generation_postResults
|
2014-09-23 16:08:20 +05:30
|
|
|
#endif
|
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
implicit none
|
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-14 16:24:45 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2014-09-23 16:08:20 +05:30
|
|
|
#ifdef multiphysicsOut
|
|
|
|
real(pReal), dimension(plasticState(material_phase(ipc,ip,el))%sizePostResults + &
|
|
|
|
damageState( material_phase(ipc,ip,el))%sizePostResults + &
|
2014-10-11 02:25:09 +05:30
|
|
|
thermalState(material_phase(ipc,ip,el))%sizePostResults + &
|
|
|
|
vacancyState(material_phase(ipc,ip,el))%sizePostResults) :: &
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults
|
|
|
|
#else
|
2014-07-04 15:16:43 +05:30
|
|
|
real(pReal), dimension(plasticState(material_phase(ipc,ip,el))%sizePostResults) :: &
|
2013-01-16 15:44:57 +05:30
|
|
|
constitutive_postResults
|
2014-09-23 16:08:20 +05:30
|
|
|
#endif
|
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) :: &
|
|
|
|
Tstar_v !< 2nd Piola Kirchhoff stress tensor (Mandel)
|
2014-10-28 08:12:25 +05:30
|
|
|
integer(pInt) :: &
|
|
|
|
startPos, endPos
|
2014-09-23 02:08:19 +05:30
|
|
|
|
2013-01-16 15:44:57 +05:30
|
|
|
constitutive_postResults = 0.0_pReal
|
|
|
|
|
2014-09-23 16:08:20 +05:30
|
|
|
startPos = 1_pInt
|
|
|
|
endPos = plasticState(material_phase(ipc,ip,el))%sizePostResults
|
2013-09-19 13:16:01 +05:30
|
|
|
select case (phase_plasticity(material_phase(ipc,ip,el)))
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_TITANMOD_ID)
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = constitutive_titanmod_postResults(ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_J2_ID)
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = constitutive_j2_postResults(Tstar_v,ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_PHENOPOWERLAW_ID)
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_phenopowerlaw_postResults(Tstar_v,ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_DISLOTWIN_ID)
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_dislotwin_postResults(Tstar_v,constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2014-08-08 16:34:40 +05:30
|
|
|
case (PLASTICITY_DISLOKMC_ID)
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_dislokmc_postResults(Tstar_v,constitutive_getTemperature(ipc,ip,el),ipc,ip,el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (PLASTICITY_NONLOCAL_ID)
|
2014-09-23 16:08:20 +05:30
|
|
|
constitutive_postResults(startPos:endPos) = &
|
2014-10-28 08:12:25 +05:30
|
|
|
constitutive_nonlocal_postResults (Tstar_v,FeArray,ip,el)
|
2013-01-16 15:44:57 +05:30
|
|
|
end select
|
2014-09-23 16:08:20 +05:30
|
|
|
|
|
|
|
#ifdef multiphysicsOut
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + damageState(material_phase(ipc,ip,el))%sizePostResults
|
|
|
|
select case (phase_damage(material_phase(ipc,ip,el)))
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_isoBrittle_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = damage_isoBrittle_postResults(ipc, ip, el)
|
|
|
|
case (LOCAL_DAMAGE_isoDuctile_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = damage_isoDuctile_postResults(ipc, ip, el)
|
2014-10-10 18:12:12 +05:30
|
|
|
case (LOCAL_DAMAGE_gurson_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = damage_gurson_postResults(ipc, ip, el)
|
2014-10-28 16:19:12 +05:30
|
|
|
case (LOCAL_DAMAGE_anisoBrittle_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = damage_anisoBrittle_postResults(ipc, ip, el)
|
2014-09-23 16:08:20 +05:30
|
|
|
end select
|
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + thermalState(material_phase(ipc,ip,el))%sizePostResults
|
|
|
|
select case (phase_thermal(material_phase(ipc,ip,el)))
|
2014-10-09 19:38:32 +05:30
|
|
|
case (LOCAL_THERMAL_ADIABATIC_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = thermal_adiabatic_postResults(ipc, ip, el)
|
2014-09-23 16:08:20 +05:30
|
|
|
end select
|
2014-10-11 02:25:09 +05:30
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + vacancyState(material_phase(ipc,ip,el))%sizePostResults
|
|
|
|
select case (phase_vacancy(material_phase(ipc,ip,el)))
|
|
|
|
case (LOCAL_VACANCY_generation_ID)
|
|
|
|
constitutive_postResults(startPos:endPos) = vacancy_generation_postResults(ipc, ip, el)
|
|
|
|
end select
|
2014-09-23 16:08:20 +05:30
|
|
|
#endif
|
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
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
end function constitutive_postResults
|
2009-03-06 15:32:36 +05:30
|
|
|
|
2012-03-14 21:46:11 +05:30
|
|
|
|
2013-02-11 16:13:45 +05:30
|
|
|
end module constitutive
|