2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Denny Tjahjanto, Max-Planck-Institut für Eisenforschung GmbH
|
2015-10-15 00:06:19 +05:30
|
|
|
!> @brief homogenization manager, organizing deformation partitioning and stress homogenization
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
module homogenization
|
2013-01-29 15:58:01 +05:30
|
|
|
use prec, only: &
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
|
|
|
tOutputData, &
|
|
|
|
#endif
|
2013-01-29 15:58:01 +05:30
|
|
|
pInt, &
|
2014-09-19 23:29:06 +05:30
|
|
|
pReal
|
2014-08-21 23:18:20 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! General variables for the homogenization at a material point
|
2009-05-07 21:57:36 +05:30
|
|
|
implicit none
|
2013-01-29 15:58:01 +05:30
|
|
|
private
|
2014-08-21 23:18:20 +05:30
|
|
|
real(pReal), dimension(:,:,:,:), allocatable, public :: &
|
2013-01-29 15:58:01 +05:30
|
|
|
materialpoint_F0, & !< def grad of IP at start of FE increment
|
|
|
|
materialpoint_F, & !< def grad of IP to be reached at end of FE increment
|
|
|
|
materialpoint_P !< first P--K stress of IP
|
|
|
|
real(pReal), dimension(:,:,:,:,:,:), allocatable, public :: &
|
|
|
|
materialpoint_dPdF !< tangent of first P--K stress at IP
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
|
|
|
type(tOutputData), dimension(:), allocatable, public :: &
|
|
|
|
homogOutput
|
|
|
|
type(tOutputData), dimension(:,:), allocatable, public :: &
|
|
|
|
crystalliteOutput, &
|
|
|
|
phaseOutput
|
|
|
|
#else
|
2013-01-29 15:58:01 +05:30
|
|
|
real(pReal), dimension(:,:,:), allocatable, public :: &
|
|
|
|
materialpoint_results !< results array of material point
|
2015-03-30 15:15:10 +05:30
|
|
|
#endif
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt), public, protected :: &
|
|
|
|
materialpoint_sizeResults, &
|
2014-09-26 16:04:36 +05:30
|
|
|
homogenization_maxSizePostResults, &
|
2015-05-28 22:32:23 +05:30
|
|
|
thermal_maxSizePostResults, &
|
|
|
|
damage_maxSizePostResults, &
|
|
|
|
vacancyflux_maxSizePostResults, &
|
|
|
|
porosity_maxSizePostResults, &
|
|
|
|
hydrogenflux_maxSizePostResults
|
2013-01-29 15:58:01 +05:30
|
|
|
|
|
|
|
real(pReal), dimension(:,:,:,:), allocatable, private :: &
|
|
|
|
materialpoint_subF0, & !< def grad of IP at beginning of homogenization increment
|
|
|
|
materialpoint_subF !< def grad of IP to be reached at end of homog inc
|
|
|
|
real(pReal), dimension(:,:), allocatable, private :: &
|
|
|
|
materialpoint_subFrac, &
|
|
|
|
materialpoint_subStep, &
|
2013-10-19 02:26:10 +05:30
|
|
|
materialpoint_subdt
|
2013-01-29 15:58:01 +05:30
|
|
|
logical, dimension(:,:), allocatable, private :: &
|
|
|
|
materialpoint_requested, &
|
|
|
|
materialpoint_converged
|
|
|
|
logical, dimension(:,:,:), allocatable, private :: &
|
|
|
|
materialpoint_doneAndHappy
|
|
|
|
|
|
|
|
public :: &
|
|
|
|
homogenization_init, &
|
|
|
|
materialpoint_stressAndItsTangent, &
|
2015-05-28 22:32:23 +05:30
|
|
|
materialpoint_postResults
|
2013-01-29 15:58:01 +05:30
|
|
|
private :: &
|
|
|
|
homogenization_partitionDeformation, &
|
|
|
|
homogenization_updateState, &
|
|
|
|
homogenization_averageStressAndItsTangent, &
|
|
|
|
homogenization_postResults
|
2012-08-10 21:28:17 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief module initialization
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-07-24 20:27:29 +05:30
|
|
|
subroutine homogenization_init
|
2012-08-10 21:28:17 +05:30
|
|
|
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
|
2013-01-29 15:58:01 +05:30
|
|
|
use math, only: &
|
|
|
|
math_I3
|
|
|
|
use debug, only: &
|
|
|
|
debug_level, &
|
|
|
|
debug_homogenization, &
|
2013-10-16 18:34:59 +05:30
|
|
|
debug_levelBasic, &
|
|
|
|
debug_e, &
|
|
|
|
debug_g
|
2013-01-29 15:58:01 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_maxNips, &
|
|
|
|
mesh_NcpElems, &
|
2015-10-15 00:06:19 +05:30
|
|
|
mesh_element, &
|
2013-01-29 15:58:01 +05:30
|
|
|
FE_Nips, &
|
|
|
|
FE_geomtype
|
2015-04-21 19:40:34 +05:30
|
|
|
#ifdef FEM
|
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_sizePostResults
|
|
|
|
#else
|
2013-01-29 15:58:01 +05:30
|
|
|
use constitutive, only: &
|
2015-05-28 22:32:23 +05:30
|
|
|
constitutive_plasticity_maxSizePostResults, &
|
|
|
|
constitutive_source_maxSizePostResults
|
2013-01-29 15:58:01 +05:30
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_maxSizePostResults
|
2015-04-21 19:40:34 +05:30
|
|
|
#endif
|
2012-08-10 21:28:17 +05:30
|
|
|
use material
|
2014-03-14 04:50:50 +05:30
|
|
|
use homogenization_none
|
2012-08-10 21:28:17 +05:30
|
|
|
use homogenization_isostrain
|
|
|
|
use homogenization_RGC
|
2015-05-28 22:32:23 +05:30
|
|
|
use thermal_isothermal
|
|
|
|
use thermal_adiabatic
|
|
|
|
use thermal_conduction
|
|
|
|
use damage_none
|
|
|
|
use damage_local
|
|
|
|
use damage_nonlocal
|
|
|
|
use vacancyflux_isoconc
|
|
|
|
use vacancyflux_isochempot
|
|
|
|
use vacancyflux_cahnhilliard
|
|
|
|
use porosity_none
|
|
|
|
use porosity_phasefield
|
|
|
|
use hydrogenflux_isoconc
|
|
|
|
use hydrogenflux_cahnhilliard
|
2014-09-26 16:04:36 +05:30
|
|
|
use IO
|
2014-10-10 01:53:06 +05:30
|
|
|
use numerics, only: &
|
|
|
|
worldrank
|
2012-05-08 20:27:06 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
implicit none
|
2013-12-12 22:39:59 +05:30
|
|
|
integer(pInt), parameter :: FILEUNIT = 200_pInt
|
2015-05-28 22:32:23 +05:30
|
|
|
integer(pInt) :: e,i,p
|
2012-08-10 21:28:17 +05:30
|
|
|
integer(pInt), dimension(:,:), pointer :: thisSize
|
2014-09-26 16:04:36 +05:30
|
|
|
integer(pInt), dimension(:) , pointer :: thisNoutput
|
2012-08-10 21:28:17 +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
|
2015-05-28 22:32:23 +05:30
|
|
|
logical :: knownHomogenization, knownThermal, knownDamage, knownVacancyflux, knownPorosity, knownHydrogenflux
|
2014-09-26 16:04:36 +05:30
|
|
|
|
2016-01-17 16:44:06 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2016-01-17 16:44:06 +05:30
|
|
|
! open material.config
|
2013-12-12 22:39:59 +05:30
|
|
|
if (.not. IO_open_jobFile_stat(FILEUNIT,material_localFileExt)) & ! no local material configuration present...
|
|
|
|
call IO_open_file(FILEUNIT,material_configFile) ! ... open material.config file
|
2016-01-17 16:44:06 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse homogenization from config file
|
2014-03-14 04:50:50 +05:30
|
|
|
if (any(homogenization_type == HOMOGENIZATION_NONE_ID)) &
|
2014-05-08 23:14:28 +05:30
|
|
|
call homogenization_none_init()
|
2014-03-14 04:50:50 +05:30
|
|
|
if (any(homogenization_type == HOMOGENIZATION_ISOSTRAIN_ID)) &
|
|
|
|
call homogenization_isostrain_init(FILEUNIT)
|
|
|
|
if (any(homogenization_type == HOMOGENIZATION_RGC_ID)) &
|
|
|
|
call homogenization_RGC_init(FILEUNIT)
|
2015-09-24 23:20:11 +05:30
|
|
|
|
2014-09-26 16:04:36 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-05-28 22:32:23 +05:30
|
|
|
! parse thermal from config file
|
2015-09-24 23:20:11 +05:30
|
|
|
call IO_checkAndRewind(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(thermal_type == THERMAL_isothermal_ID)) &
|
2015-07-24 20:23:50 +05:30
|
|
|
call thermal_isothermal_init()
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(thermal_type == THERMAL_adiabatic_ID)) &
|
2015-07-24 20:23:50 +05:30
|
|
|
call thermal_adiabatic_init(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(thermal_type == THERMAL_conduction_ID)) &
|
2015-07-24 20:23:50 +05:30
|
|
|
call thermal_conduction_init(FILEUNIT)
|
2015-09-24 23:20:11 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse damage from config file
|
2015-09-24 23:20:11 +05:30
|
|
|
call IO_checkAndRewind(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(damage_type == DAMAGE_none_ID)) &
|
|
|
|
call damage_none_init()
|
|
|
|
if (any(damage_type == DAMAGE_local_ID)) &
|
|
|
|
call damage_local_init(FILEUNIT)
|
|
|
|
if (any(damage_type == DAMAGE_nonlocal_ID)) &
|
|
|
|
call damage_nonlocal_init(FILEUNIT)
|
2015-09-24 23:20:11 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse vacancy transport from config file
|
2015-09-24 23:20:11 +05:30
|
|
|
call IO_checkAndRewind(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(vacancyflux_type == VACANCYFLUX_isoconc_ID)) &
|
|
|
|
call vacancyflux_isoconc_init()
|
|
|
|
if (any(vacancyflux_type == VACANCYFLUX_isochempot_ID)) &
|
|
|
|
call vacancyflux_isochempot_init(FILEUNIT)
|
|
|
|
if (any(vacancyflux_type == VACANCYFLUX_cahnhilliard_ID)) &
|
|
|
|
call vacancyflux_cahnhilliard_init(FILEUNIT)
|
2015-09-24 23:20:11 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse porosity from config file
|
2015-09-24 23:20:11 +05:30
|
|
|
call IO_checkAndRewind(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(porosity_type == POROSITY_none_ID)) &
|
|
|
|
call porosity_none_init()
|
|
|
|
if (any(porosity_type == POROSITY_phasefield_ID)) &
|
|
|
|
call porosity_phasefield_init(FILEUNIT)
|
2015-09-24 23:20:11 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! parse hydrogen transport from config file
|
2015-09-24 23:20:11 +05:30
|
|
|
call IO_checkAndRewind(FILEUNIT)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (any(hydrogenflux_type == HYDROGENFLUX_isoconc_ID)) &
|
|
|
|
call hydrogenflux_isoconc_init()
|
|
|
|
if (any(hydrogenflux_type == HYDROGENFLUX_cahnhilliard_ID)) &
|
|
|
|
call hydrogenflux_cahnhilliard_init(FILEUNIT)
|
2014-09-26 16:04:36 +05:30
|
|
|
close(FILEUNIT)
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! write description file for homogenization output
|
2015-08-18 21:37:01 +05:30
|
|
|
mainProcess2: if (worldrank == 0) then
|
2015-03-25 21:32:30 +05:30
|
|
|
call IO_write_jobFile(FILEUNIT,'outputHomogenization')
|
|
|
|
do p = 1,material_Nhomogenization
|
2015-04-21 19:40:34 +05:30
|
|
|
if (any(material_homog == p)) then
|
|
|
|
i = homogenization_typeInstance(p) ! which instance of this homogenization type
|
|
|
|
knownHomogenization = .true. ! assume valid
|
|
|
|
select case(homogenization_type(p)) ! split per homogenization type
|
|
|
|
case (HOMOGENIZATION_NONE_ID)
|
|
|
|
outputName = HOMOGENIZATION_NONE_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (HOMOGENIZATION_ISOSTRAIN_ID)
|
|
|
|
outputName = HOMOGENIZATION_ISOSTRAIN_label
|
|
|
|
thisNoutput => homogenization_isostrain_Noutput
|
|
|
|
thisOutput => homogenization_isostrain_output
|
|
|
|
thisSize => homogenization_isostrain_sizePostResult
|
|
|
|
case (HOMOGENIZATION_RGC_ID)
|
|
|
|
outputName = HOMOGENIZATION_RGC_label
|
|
|
|
thisNoutput => homogenization_RGC_Noutput
|
|
|
|
thisOutput => homogenization_RGC_output
|
|
|
|
thisSize => homogenization_RGC_sizePostResult
|
|
|
|
case default
|
|
|
|
knownHomogenization = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-04-21 19:40:34 +05:30
|
|
|
write(FILEUNIT,'(/,a,/)') '['//trim(homogenization_name(p))//']'
|
|
|
|
if (knownHomogenization) then
|
|
|
|
write(FILEUNIT,'(a)') '(type)'//char(9)//trim(outputName)
|
|
|
|
write(FILEUNIT,'(a,i4)') '(ngrains)'//char(9),homogenization_Ngrains(p)
|
|
|
|
if (homogenization_type(p) /= HOMOGENIZATION_NONE_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-05-28 22:32:23 +05:30
|
|
|
i = thermal_typeInstance(p) ! which instance of this thermal type
|
|
|
|
knownThermal = .true. ! assume valid
|
|
|
|
select case(thermal_type(p)) ! split per thermal type
|
|
|
|
case (THERMAL_isothermal_ID)
|
|
|
|
outputName = THERMAL_isothermal_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (THERMAL_adiabatic_ID)
|
|
|
|
outputName = THERMAL_adiabatic_label
|
|
|
|
thisNoutput => thermal_adiabatic_Noutput
|
|
|
|
thisOutput => thermal_adiabatic_output
|
|
|
|
thisSize => thermal_adiabatic_sizePostResult
|
|
|
|
case (THERMAL_conduction_ID)
|
|
|
|
outputName = THERMAL_conduction_label
|
|
|
|
thisNoutput => thermal_conduction_Noutput
|
|
|
|
thisOutput => thermal_conduction_output
|
|
|
|
thisSize => thermal_conduction_sizePostResult
|
|
|
|
case default
|
|
|
|
knownThermal = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-05-28 22:32:23 +05:30
|
|
|
if (knownThermal) then
|
|
|
|
write(FILEUNIT,'(a)') '(thermal)'//char(9)//trim(outputName)
|
|
|
|
if (thermal_type(p) /= THERMAL_isothermal_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-05-28 22:32:23 +05:30
|
|
|
i = damage_typeInstance(p) ! which instance of this damage type
|
|
|
|
knownDamage = .true. ! assume valid
|
|
|
|
select case(damage_type(p)) ! split per damage type
|
|
|
|
case (DAMAGE_none_ID)
|
|
|
|
outputName = DAMAGE_none_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (DAMAGE_local_ID)
|
|
|
|
outputName = DAMAGE_local_label
|
|
|
|
thisNoutput => damage_local_Noutput
|
|
|
|
thisOutput => damage_local_output
|
|
|
|
thisSize => damage_local_sizePostResult
|
|
|
|
case (DAMAGE_nonlocal_ID)
|
|
|
|
outputName = DAMAGE_nonlocal_label
|
|
|
|
thisNoutput => damage_nonlocal_Noutput
|
|
|
|
thisOutput => damage_nonlocal_output
|
|
|
|
thisSize => damage_nonlocal_sizePostResult
|
|
|
|
case default
|
|
|
|
knownDamage = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-05-28 22:32:23 +05:30
|
|
|
if (knownDamage) then
|
|
|
|
write(FILEUNIT,'(a)') '(damage)'//char(9)//trim(outputName)
|
|
|
|
if (damage_type(p) /= DAMAGE_none_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-05-28 22:32:23 +05:30
|
|
|
i = vacancyflux_typeInstance(p) ! which instance of this vacancy flux type
|
|
|
|
knownVacancyflux = .true. ! assume valid
|
|
|
|
select case(vacancyflux_type(p)) ! split per vacancy flux type
|
|
|
|
case (VACANCYFLUX_isoconc_ID)
|
|
|
|
outputName = VACANCYFLUX_isoconc_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (VACANCYFLUX_isochempot_ID)
|
|
|
|
outputName = VACANCYFLUX_isochempot_label
|
|
|
|
thisNoutput => vacancyflux_isochempot_Noutput
|
|
|
|
thisOutput => vacancyflux_isochempot_output
|
|
|
|
thisSize => vacancyflux_isochempot_sizePostResult
|
|
|
|
case (VACANCYFLUX_cahnhilliard_ID)
|
|
|
|
outputName = VACANCYFLUX_cahnhilliard_label
|
|
|
|
thisNoutput => vacancyflux_cahnhilliard_Noutput
|
|
|
|
thisOutput => vacancyflux_cahnhilliard_output
|
|
|
|
thisSize => vacancyflux_cahnhilliard_sizePostResult
|
|
|
|
case default
|
|
|
|
knownVacancyflux = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-05-28 22:32:23 +05:30
|
|
|
if (knownVacancyflux) then
|
|
|
|
write(FILEUNIT,'(a)') '(vacancyflux)'//char(9)//trim(outputName)
|
|
|
|
if (vacancyflux_type(p) /= VACANCYFLUX_isoconc_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-05-28 22:32:23 +05:30
|
|
|
i = porosity_typeInstance(p) ! which instance of this porosity type
|
|
|
|
knownPorosity = .true. ! assume valid
|
|
|
|
select case(porosity_type(p)) ! split per porosity type
|
|
|
|
case (POROSITY_none_ID)
|
|
|
|
outputName = POROSITY_none_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (POROSITY_phasefield_ID)
|
|
|
|
outputName = POROSITY_phasefield_label
|
|
|
|
thisNoutput => porosity_phasefield_Noutput
|
|
|
|
thisOutput => porosity_phasefield_output
|
|
|
|
thisSize => porosity_phasefield_sizePostResult
|
|
|
|
case default
|
|
|
|
knownPorosity = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-05-28 22:32:23 +05:30
|
|
|
if (knownPorosity) then
|
|
|
|
write(FILEUNIT,'(a)') '(porosity)'//char(9)//trim(outputName)
|
|
|
|
if (porosity_type(p) /= POROSITY_none_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-05-28 22:32:23 +05:30
|
|
|
i = hydrogenflux_typeInstance(p) ! which instance of this hydrogen flux type
|
|
|
|
knownHydrogenflux = .true. ! assume valid
|
|
|
|
select case(hydrogenflux_type(p)) ! split per hydrogen flux type
|
|
|
|
case (HYDROGENFLUX_isoconc_ID)
|
|
|
|
outputName = HYDROGENFLUX_isoconc_label
|
|
|
|
thisNoutput => null()
|
|
|
|
thisOutput => null()
|
|
|
|
thisSize => null()
|
|
|
|
case (HYDROGENFLUX_cahnhilliard_ID)
|
|
|
|
outputName = HYDROGENFLUX_cahnhilliard_label
|
|
|
|
thisNoutput => hydrogenflux_cahnhilliard_Noutput
|
|
|
|
thisOutput => hydrogenflux_cahnhilliard_output
|
|
|
|
thisSize => hydrogenflux_cahnhilliard_sizePostResult
|
|
|
|
case default
|
|
|
|
knownHydrogenflux = .false.
|
2015-10-15 00:06:19 +05:30
|
|
|
end select
|
2015-05-28 22:32:23 +05:30
|
|
|
if (knownHydrogenflux) then
|
|
|
|
write(FILEUNIT,'(a)') '(hydrogenflux)'//char(9)//trim(outputName)
|
|
|
|
if (hydrogenflux_type(p) /= HYDROGENFLUX_isoconc_ID) then
|
|
|
|
do e = 1,thisNoutput(i)
|
|
|
|
write(FILEUNIT,'(a,i4)') trim(thisOutput(e,i))//char(9),thisSize(e,i)
|
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
|
|
|
endif
|
2015-04-21 19:40:34 +05:30
|
|
|
endif
|
2014-09-26 16:04:36 +05:30
|
|
|
enddo
|
2015-03-25 21:32:30 +05:30
|
|
|
close(FILEUNIT)
|
2015-08-18 21:37:01 +05:30
|
|
|
endif mainProcess2
|
2014-09-22 23:45:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! allocate and initialize global variables
|
2013-12-16 17:28:03 +05:30
|
|
|
allocate(materialpoint_dPdF(3,3,3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_F0(3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
2013-01-29 15:58:01 +05:30
|
|
|
materialpoint_F0 = spread(spread(math_I3,3,mesh_maxNips),4,mesh_NcpElems) ! initialize to identity
|
2013-12-16 17:28:03 +05:30
|
|
|
allocate(materialpoint_F(3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
materialpoint_F = materialpoint_F0 ! initialize to identity
|
|
|
|
allocate(materialpoint_subF0(3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_subF(3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_P(3,3,mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_subFrac(mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_subStep(mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_subdt(mesh_maxNips,mesh_NcpElems), source=0.0_pReal)
|
|
|
|
allocate(materialpoint_requested(mesh_maxNips,mesh_NcpElems), source=.false.)
|
|
|
|
allocate(materialpoint_converged(mesh_maxNips,mesh_NcpElems), source=.true.)
|
|
|
|
allocate(materialpoint_doneAndHappy(2,mesh_maxNips,mesh_NcpElems), source=.true.)
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-16 18:34:59 +05:30
|
|
|
! allocate and initialize global state and postresutls variables
|
2014-09-26 16:02:56 +05:30
|
|
|
homogenization_maxSizePostResults = 0_pInt
|
2015-05-28 22:32:23 +05:30
|
|
|
thermal_maxSizePostResults = 0_pInt
|
|
|
|
damage_maxSizePostResults = 0_pInt
|
|
|
|
vacancyflux_maxSizePostResults = 0_pInt
|
|
|
|
porosity_maxSizePostResults = 0_pInt
|
|
|
|
hydrogenflux_maxSizePostResults = 0_pInt
|
2014-09-26 16:02:56 +05:30
|
|
|
do p = 1,material_Nhomogenization
|
2015-05-28 22:32:23 +05:30
|
|
|
homogenization_maxSizePostResults = max(homogenization_maxSizePostResults,homogState (p)%sizePostResults)
|
|
|
|
thermal_maxSizePostResults = max(thermal_maxSizePostResults, thermalState (p)%sizePostResults)
|
|
|
|
damage_maxSizePostResults = max(damage_maxSizePostResults ,damageState (p)%sizePostResults)
|
|
|
|
vacancyflux_maxSizePostResults = max(vacancyflux_maxSizePostResults ,vacancyfluxState (p)%sizePostResults)
|
|
|
|
porosity_maxSizePostResults = max(porosity_maxSizePostResults ,porosityState (p)%sizePostResults)
|
|
|
|
hydrogenflux_maxSizePostResults = max(hydrogenflux_maxSizePostResults ,hydrogenfluxState(p)%sizePostResults)
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo
|
2015-03-30 15:15:10 +05:30
|
|
|
|
|
|
|
#ifdef FEM
|
|
|
|
allocate(homogOutput (material_Nhomogenization ))
|
|
|
|
allocate(crystalliteOutput(material_Ncrystallite, homogenization_maxNgrains))
|
|
|
|
allocate(phaseOutput (material_Nphase, homogenization_maxNgrains))
|
2015-10-15 00:06:19 +05:30
|
|
|
do p = 1, material_Nhomogenization
|
2015-05-28 22:32:23 +05:30
|
|
|
homogOutput(p)%sizeResults = homogState (p)%sizePostResults + &
|
|
|
|
thermalState (p)%sizePostResults + &
|
|
|
|
damageState (p)%sizePostResults + &
|
|
|
|
vacancyfluxState (p)%sizePostResults + &
|
|
|
|
porosityState (p)%sizePostResults + &
|
|
|
|
hydrogenfluxState(p)%sizePostResults
|
2015-03-30 15:15:10 +05:30
|
|
|
homogOutput(p)%sizeIpCells = count(material_homog==p)
|
|
|
|
allocate(homogOutput(p)%output(homogOutput(p)%sizeResults,homogOutput(p)%sizeIpCells))
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo
|
2015-03-30 15:15:10 +05:30
|
|
|
do p = 1, material_Ncrystallite; do e = 1, homogenization_maxNgrains
|
|
|
|
crystalliteOutput(p,e)%sizeResults = crystallite_sizePostResults(p)
|
|
|
|
crystalliteOutput(p,e)%sizeIpCells = count(microstructure_crystallite(mesh_element(4,:)) == p .and. &
|
|
|
|
homogenization_Ngrains (mesh_element(3,:)) >= e)*mesh_maxNips
|
|
|
|
allocate(crystalliteOutput(p,e)%output(crystalliteOutput(p,e)%sizeResults,crystalliteOutput(p,e)%sizeIpCells))
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo; enddo
|
2015-03-30 15:15:10 +05:30
|
|
|
do p = 1, material_Nphase; do e = 1, homogenization_maxNgrains
|
2015-05-28 22:32:23 +05:30
|
|
|
phaseOutput(p,e)%sizeResults = plasticState (p)%sizePostResults + &
|
|
|
|
sum(sourceState (p)%p(:)%sizePostResults)
|
2015-03-30 15:15:10 +05:30
|
|
|
phaseOutput(p,e)%sizeIpCells = count(material_phase(e,:,:) == p)
|
|
|
|
allocate(phaseOutput(p,e)%output(phaseOutput(p,e)%sizeResults,phaseOutput(p,e)%sizeIpCells))
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo; enddo
|
2015-03-30 15:15:10 +05:30
|
|
|
#else
|
2013-01-29 15:58:01 +05:30
|
|
|
materialpoint_sizeResults = 1 & ! grain count
|
|
|
|
+ 1 + homogenization_maxSizePostResults & ! homogSize & homogResult
|
2015-05-28 22:32:23 +05:30
|
|
|
+ thermal_maxSizePostResults &
|
|
|
|
+ damage_maxSizePostResults &
|
|
|
|
+ vacancyflux_maxSizePostResults &
|
|
|
|
+ porosity_maxSizePostResults &
|
|
|
|
+ hydrogenflux_maxSizePostResults &
|
2013-01-29 15:58:01 +05:30
|
|
|
+ homogenization_maxNgrains * (1 + crystallite_maxSizePostResults & ! crystallite size & crystallite results
|
2015-05-28 22:32:23 +05:30
|
|
|
+ 1 + constitutive_plasticity_maxSizePostResults & ! constitutive size & constitutive results
|
2015-10-15 00:06:19 +05:30
|
|
|
+ constitutive_source_maxSizePostResults)
|
2012-08-10 21:28:17 +05:30
|
|
|
allocate(materialpoint_results(materialpoint_sizeResults,mesh_maxNips,mesh_NcpElems))
|
2015-03-30 15:15:10 +05:30
|
|
|
#endif
|
2015-10-15 00:06:19 +05:30
|
|
|
|
|
|
|
mainProcess: if (worldrank == 0) then
|
2014-10-10 21:51:10 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- homogenization init -+>>>'
|
|
|
|
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
|
2012-02-01 00:48:55 +05:30
|
|
|
#include "compilation_info.f90"
|
2014-10-10 21:51:10 +05:30
|
|
|
endif mainProcess
|
2014-10-10 01:53:06 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt) then
|
2014-09-19 23:29:06 +05:30
|
|
|
#ifdef TODO
|
2013-10-16 18:34:59 +05:30
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'homogenization_state0: ', shape(homogenization_state0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'homogenization_subState0: ', shape(homogenization_subState0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'homogenization_state: ', shape(homogenization_state)
|
2014-08-21 23:18:20 +05:30
|
|
|
#endif
|
2013-10-16 18:34:59 +05:30
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_dPdF: ', shape(materialpoint_dPdF)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_F0: ', shape(materialpoint_F0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_F: ', shape(materialpoint_F)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subF0: ', shape(materialpoint_subF0)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subF: ', shape(materialpoint_subF)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_P: ', shape(materialpoint_P)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subFrac: ', shape(materialpoint_subFrac)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subStep: ', shape(materialpoint_subStep)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subdt: ', shape(materialpoint_subdt)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_requested: ', shape(materialpoint_requested)
|
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'materialpoint_converged: ', shape(materialpoint_converged)
|
|
|
|
write(6,'(a32,1x,7(i8,1x),/)') 'materialpoint_doneAndHappy: ', shape(materialpoint_doneAndHappy)
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifndef FEM
|
2013-10-16 18:34:59 +05:30
|
|
|
write(6,'(a32,1x,7(i8,1x),/)') 'materialpoint_results: ', shape(materialpoint_results)
|
2015-03-30 15:15:10 +05:30
|
|
|
#endif
|
2013-10-16 18:34:59 +05:30
|
|
|
write(6,'(a32,1x,7(i8,1x))') 'maxSizePostResults: ', homogenization_maxSizePostResults
|
2013-01-29 15:58:01 +05:30
|
|
|
endif
|
|
|
|
flush(6)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-10-16 18:34:59 +05:30
|
|
|
if (debug_g < 1 .or. debug_g > homogenization_Ngrains(mesh_element(3,debug_e))) &
|
|
|
|
call IO_error(602_pInt,ext_msg='component (grain)')
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end subroutine homogenization_init
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief parallelized calculation of stress and corresponding tangent at material points
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine materialpoint_stressAndItsTangent(updateJaco,dt)
|
2013-01-29 15:58:01 +05:30
|
|
|
use numerics, only: &
|
|
|
|
subStepMinHomog, &
|
|
|
|
subStepSizeHomog, &
|
|
|
|
stepIncreaseHomog, &
|
|
|
|
nHomog, &
|
|
|
|
nMPstate
|
2015-10-15 00:06:19 +05:30
|
|
|
use math, only: &
|
2013-01-29 15:58:01 +05:30
|
|
|
math_transpose33
|
|
|
|
use FEsolving, only: &
|
|
|
|
FEsolving_execElem, &
|
|
|
|
FEsolving_execIP, &
|
|
|
|
terminallyIll
|
|
|
|
use mesh, only: &
|
2015-04-21 19:40:34 +05:30
|
|
|
mesh_element
|
2013-01-29 15:58:01 +05:30
|
|
|
use material, only: &
|
2014-05-12 06:14:44 +05:30
|
|
|
plasticState, &
|
2015-05-28 22:32:23 +05:30
|
|
|
sourceState, &
|
2014-08-21 23:18:20 +05:30
|
|
|
homogState, &
|
2015-05-28 22:32:23 +05:30
|
|
|
thermalState, &
|
|
|
|
damageState, &
|
|
|
|
vacancyfluxState, &
|
|
|
|
porosityState, &
|
|
|
|
hydrogenfluxState, &
|
|
|
|
phase_Nsources, &
|
2015-10-15 00:06:19 +05:30
|
|
|
mappingHomogenization, &
|
2016-01-15 05:49:44 +05:30
|
|
|
phaseAt, phasememberAt, &
|
2013-01-29 15:58:01 +05:30
|
|
|
homogenization_Ngrains
|
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_F0, &
|
|
|
|
crystallite_Fp0, &
|
|
|
|
crystallite_Fp, &
|
2015-03-06 18:39:00 +05:30
|
|
|
crystallite_Fi0, &
|
|
|
|
crystallite_Fi, &
|
2013-01-29 15:58:01 +05:30
|
|
|
crystallite_Lp0, &
|
|
|
|
crystallite_Lp, &
|
2015-03-06 18:39:00 +05:30
|
|
|
crystallite_Li0, &
|
|
|
|
crystallite_Li, &
|
2013-01-29 15:58:01 +05:30
|
|
|
crystallite_dPdF, &
|
|
|
|
crystallite_dPdF0, &
|
|
|
|
crystallite_Tstar0_v, &
|
|
|
|
crystallite_Tstar_v, &
|
|
|
|
crystallite_partionedF0, &
|
|
|
|
crystallite_partionedF, &
|
|
|
|
crystallite_partionedFp0, &
|
|
|
|
crystallite_partionedLp0, &
|
2015-03-06 18:39:00 +05:30
|
|
|
crystallite_partionedFi0, &
|
|
|
|
crystallite_partionedLi0, &
|
2013-01-29 15:58:01 +05:30
|
|
|
crystallite_partioneddPdF0, &
|
|
|
|
crystallite_partionedTstar0_v, &
|
|
|
|
crystallite_dt, &
|
|
|
|
crystallite_requested, &
|
|
|
|
crystallite_converged, &
|
|
|
|
crystallite_stressAndItsTangent, &
|
|
|
|
crystallite_orientations
|
|
|
|
use debug, only: &
|
|
|
|
debug_level, &
|
|
|
|
debug_homogenization, &
|
|
|
|
debug_levelBasic, &
|
|
|
|
debug_levelSelective, &
|
|
|
|
debug_e, &
|
|
|
|
debug_i, &
|
|
|
|
debug_MaterialpointLoopDistribution, &
|
|
|
|
debug_MaterialpointStateLoopDistribution
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2009-05-07 21:57:36 +05:30
|
|
|
implicit none
|
2012-08-10 21:28:17 +05:30
|
|
|
real(pReal), intent(in) :: dt !< time increment
|
|
|
|
logical, intent(in) :: updateJaco !< initiating Jacobian update
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt) :: &
|
|
|
|
NiterationHomog, &
|
|
|
|
NiterationMPstate, &
|
|
|
|
g, & !< grain number
|
|
|
|
i, & !< integration point number
|
|
|
|
e, & !< element number
|
2015-05-28 22:32:23 +05:30
|
|
|
mySource, &
|
2013-01-29 15:58:01 +05:30
|
|
|
myNgrains
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! initialize to starting condition
|
2013-10-16 18:34:59 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt) then
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP CRITICAL (write2out)
|
2013-10-16 18:34:59 +05:30
|
|
|
write(6,'(/a,i5,1x,i2)') '<< HOMOG >> Material Point start at el ip ', debug_e, debug_i
|
|
|
|
|
2012-10-02 18:23:25 +05:30
|
|
|
write(6,'(a,/,3(12x,3(f14.9,1x)/))') '<< HOMOG >> F0', &
|
2012-08-10 21:28:17 +05:30
|
|
|
math_transpose33(materialpoint_F0(1:3,1:3,debug_i,debug_e))
|
2012-10-02 18:23:25 +05:30
|
|
|
write(6,'(a,/,3(12x,3(f14.9,1x)/))') '<< HOMOG >> F', &
|
2012-08-10 21:28:17 +05:30
|
|
|
math_transpose33(materialpoint_F(1:3,1:3,debug_i,debug_e))
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP END CRITICAL (write2out)
|
2010-03-19 19:44:08 +05:30
|
|
|
endif
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-11-07 21:13:29 +05:30
|
|
|
! initialize restoration points of ...
|
|
|
|
do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
2009-05-07 21:57:36 +05:30
|
|
|
myNgrains = homogenization_Ngrains(mesh_element(3,e))
|
2015-05-28 22:32:23 +05:30
|
|
|
do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e); do g = 1,myNgrains
|
2014-09-03 01:46:33 +05:30
|
|
|
|
2016-01-15 05:49:44 +05:30
|
|
|
plasticState (phaseAt(g,i,e))%partionedState0(:,phasememberAt(g,i,e)) = &
|
|
|
|
plasticState (phaseAt(g,i,e))%state0( :,phasememberAt(g,i,e))
|
|
|
|
do mySource = 1_pInt, phase_Nsources(phaseAt(g,i,e))
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%partionedState0(:,phasememberAt(g,i,e)) = &
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%state0( :,phasememberAt(g,i,e))
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo
|
2014-09-03 01:46:33 +05:30
|
|
|
|
2012-11-07 21:13:29 +05:30
|
|
|
crystallite_partionedFp0(1:3,1:3,g,i,e) = crystallite_Fp0(1:3,1:3,g,i,e) ! ...plastic def grads
|
2016-01-09 20:33:18 +05:30
|
|
|
crystallite_partionedLp0(1:3,1:3,g,i,e) = crystallite_Lp0(1:3,1:3,g,i,e) ! ...plastic velocity grads
|
2015-03-06 18:39:00 +05:30
|
|
|
crystallite_partionedFi0(1:3,1:3,g,i,e) = crystallite_Fi0(1:3,1:3,g,i,e) ! ...intermediate def grads
|
|
|
|
crystallite_partionedLi0(1:3,1:3,g,i,e) = crystallite_Li0(1:3,1:3,g,i,e) ! ...intermediate velocity grads
|
2012-11-07 21:13:29 +05:30
|
|
|
crystallite_partioneddPdF0(1:3,1:3,1:3,1:3,g,i,e) = crystallite_dPdF0(1:3,1:3,1:3,1:3,g,i,e) ! ...stiffness
|
|
|
|
crystallite_partionedF0(1:3,1:3,g,i,e) = crystallite_F0(1:3,1:3,g,i,e) ! ...def grads
|
|
|
|
crystallite_partionedTstar0_v(1:6,g,i,e) = crystallite_Tstar0_v(1:6,g,i,e) ! ...2nd PK stress
|
2014-05-22 17:37:50 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
enddo; enddo
|
2012-11-07 21:13:29 +05:30
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e))
|
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
|
|
|
materialpoint_subF0(1:3,1:3,i,e) = materialpoint_F0(1:3,1:3,i,e) ! ...def grad
|
2009-05-07 21:57:36 +05:30
|
|
|
materialpoint_subFrac(i,e) = 0.0_pReal
|
2009-11-10 19:06:27 +05:30
|
|
|
materialpoint_subStep(i,e) = 1.0_pReal/subStepSizeHomog ! <<added to adopt flexibility in cutback size>>
|
2009-06-16 14:33:30 +05:30
|
|
|
materialpoint_converged(i,e) = .false. ! pretend failed step of twice the required size
|
|
|
|
materialpoint_requested(i,e) = .true. ! everybody requires calculation
|
2012-11-07 21:13:29 +05:30
|
|
|
endforall
|
2014-09-19 23:29:06 +05:30
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
2015-05-28 22:32:23 +05:30
|
|
|
homogState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal homogenization state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal thermal state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal damage state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal vacancy transport state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal porosity state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%State0( :,mappingHomogenization(1,i,e)) ! ...internal hydrogen transport state
|
2009-05-07 21:57:36 +05:30
|
|
|
enddo
|
2009-08-11 22:01:57 +05:30
|
|
|
NiterationHomog = 0_pInt
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
cutBackLooping: do while (.not. terminallyIll .and. &
|
|
|
|
any(materialpoint_subStep(:,FEsolving_execELem(1):FEsolving_execElem(2)) > subStepMinHomog))
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP PARALLEL DO PRIVATE(myNgrains)
|
2013-01-29 15:58:01 +05:30
|
|
|
elementLooping1: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
2009-05-07 21:57:36 +05:30
|
|
|
myNgrains = homogenization_Ngrains(mesh_element(3,e))
|
2013-01-29 15:58:01 +05:30
|
|
|
IpLooping1: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
converged: if ( materialpoint_converged(i,e) ) then
|
2011-03-29 12:57:19 +05:30
|
|
|
#ifndef _OPENMP
|
2012-07-05 15:24:50 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt &
|
2015-10-15 00:06:19 +05:30
|
|
|
.and. ((e == debug_e .and. i == debug_i) &
|
2012-07-05 15:24:50 +05:30
|
|
|
.or. .not. iand(debug_level(debug_homogenization),debug_levelSelective) /= 0_pInt)) then
|
2012-11-21 22:28:14 +05:30
|
|
|
write(6,'(a,1x,f12.8,1x,a,1x,f12.8,1x,a,i8,1x,i2/)') '<< HOMOG >> winding forward from', &
|
2011-03-29 12:57:19 +05:30
|
|
|
materialpoint_subFrac(i,e), 'to current materialpoint_subFrac', &
|
2012-11-21 22:28:14 +05:30
|
|
|
materialpoint_subFrac(i,e)+materialpoint_subStep(i,e),'in materialpoint_stressAndItsTangent at el ip',e,i
|
2009-08-24 13:46:01 +05:30
|
|
|
endif
|
2011-03-29 12:57:19 +05:30
|
|
|
#endif
|
2015-10-15 00:06:19 +05:30
|
|
|
|
|
|
|
!---------------------------------------------------------------------------------------------------
|
2013-01-29 15:58:01 +05:30
|
|
|
! calculate new subStep and new subFrac
|
2009-08-27 17:40:06 +05:30
|
|
|
materialpoint_subFrac(i,e) = materialpoint_subFrac(i,e) + materialpoint_subStep(i,e)
|
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
|
|
|
!$OMP FLUSH(materialpoint_subFrac)
|
2009-11-10 19:06:27 +05:30
|
|
|
materialpoint_subStep(i,e) = min(1.0_pReal-materialpoint_subFrac(i,e), &
|
2015-05-28 22:32:23 +05:30
|
|
|
stepIncreaseHomog*materialpoint_subStep(i,e)) ! introduce flexibility for step increase/acceleration
|
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
|
|
|
!$OMP FLUSH(materialpoint_subStep)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
steppingNeeded: if (materialpoint_subStep(i,e) > subStepMinHomog) then
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2009-06-16 14:33:30 +05:30
|
|
|
! wind forward grain starting point of...
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedF0(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedF(1:3,1:3,1:myNgrains,i,e) ! ...def grads
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedFp0(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_Fp(1:3,1:3,1:myNgrains,i,e) ! ...plastic def grads
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedLp0(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_Lp(1:3,1:3,1:myNgrains,i,e) ! ...plastic velocity grads
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedFi0(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_Fi(1:3,1:3,1:myNgrains,i,e) ! ...intermediate def grads
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedLi0(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_Li(1:3,1:3,1:myNgrains,i,e) ! ...intermediate velocity grads
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partioneddPdF0(1:3,1:3,1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_dPdF(1:3,1:3,1:3,1:3,1:myNgrains,i,e) ! ...stiffness
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedTstar0_v(1:6,1:myNgrains,i,e) = &
|
|
|
|
crystallite_Tstar_v(1:6,1:myNgrains,i,e) ! ...2nd PK stress
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
do g = 1,myNgrains
|
2016-01-15 05:49:44 +05:30
|
|
|
plasticState (phaseAt(g,i,e))%partionedState0(:,phasememberAt(g,i,e)) = &
|
|
|
|
plasticState (phaseAt(g,i,e))%state( :,phasememberAt(g,i,e))
|
|
|
|
do mySource = 1_pInt, phase_Nsources(phaseAt(g,i,e))
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%partionedState0(:,phasememberAt(g,i,e)) = &
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%state( :,phasememberAt(g,i,e))
|
2015-05-28 22:32:23 +05:30
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) ! ...internal homogenization state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) ! ...internal thermal state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) ! ...internal damage state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e))! ...internal vacancy transport state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e))! ...internal porosity state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) = &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e))! ...internal hydrogen transport state
|
|
|
|
materialpoint_subF0(1:3,1:3,i,e) = materialpoint_subF(1:3,1:3,i,e) ! ...def grad
|
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
|
|
|
!$OMP FLUSH(materialpoint_subF0)
|
2015-05-28 22:32:23 +05:30
|
|
|
elseif (materialpoint_requested(i,e)) then steppingNeeded ! already at final time (??)
|
2012-07-05 15:24:50 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt) then
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP CRITICAL (distributionHomog)
|
|
|
|
debug_MaterialpointLoopDistribution(min(nHomog+1,NiterationHomog)) = &
|
2014-09-03 01:46:33 +05:30
|
|
|
debug_MaterialpointLoopDistribution(min(nHomog+1,NiterationHomog)) + 1
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP END CRITICAL (distributionHomog)
|
|
|
|
endif
|
2013-01-29 15:58:01 +05:30
|
|
|
endif steppingNeeded
|
|
|
|
|
|
|
|
else converged
|
2015-05-28 22:32:23 +05:30
|
|
|
if ( (myNgrains == 1_pInt .and. materialpoint_subStep(i,e) <= 1.0 ) .or. & ! single grain already tried internal subStepping in crystallite
|
|
|
|
subStepSizeHomog * materialpoint_subStep(i,e) <= subStepMinHomog ) then ! would require too small subStep
|
|
|
|
! cutback makes no sense
|
2012-12-16 16:24:13 +05:30
|
|
|
!$OMP FLUSH(terminallyIll)
|
2015-05-28 22:32:23 +05:30
|
|
|
if (.not. terminallyIll) then ! so first signals terminally ill...
|
2012-10-18 19:18:06 +05:30
|
|
|
!$OMP CRITICAL (write2out)
|
|
|
|
write(6,*) 'Integration point ', i,' at element ', e, ' terminally ill'
|
|
|
|
!$OMP END CRITICAL (write2out)
|
|
|
|
endif
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP CRITICAL (setTerminallyIll)
|
2015-05-28 22:32:23 +05:30
|
|
|
terminallyIll = .true. ! ...and kills all others
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END CRITICAL (setTerminallyIll)
|
2015-05-28 22:32:23 +05:30
|
|
|
else ! cutback makes sense
|
|
|
|
materialpoint_subStep(i,e) = subStepSizeHomog * materialpoint_subStep(i,e) ! crystallite had severe trouble, so do a significant cutback
|
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
|
|
|
!$OMP FLUSH(materialpoint_subStep)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2011-03-29 12:57:19 +05:30
|
|
|
#ifndef _OPENMP
|
2012-07-05 15:24:50 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt &
|
2012-03-09 01:55:28 +05:30
|
|
|
.and. ((e == debug_e .and. i == debug_i) &
|
2012-07-05 15:24:50 +05:30
|
|
|
.or. .not. iand(debug_level(debug_homogenization), debug_levelSelective) /= 0_pInt)) then
|
2012-11-21 22:28:14 +05:30
|
|
|
write(6,'(a,1x,f12.8,a,i8,1x,i2/)') &
|
2011-08-02 18:06:08 +05:30
|
|
|
'<< HOMOG >> cutback step in materialpoint_stressAndItsTangent with new materialpoint_subStep:',&
|
2015-10-15 00:06:19 +05:30
|
|
|
materialpoint_subStep(i,e),' at el ip',e,i
|
2010-09-02 02:34:02 +05:30
|
|
|
endif
|
2011-03-29 12:57:19 +05:30
|
|
|
#endif
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! restore...
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_Fp(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedFp0(1:3,1:3,1:myNgrains,i,e) ! ...plastic def grads
|
|
|
|
crystallite_Lp(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedLp0(1:3,1:3,1:myNgrains,i,e) ! ...plastic velocity grads
|
|
|
|
crystallite_Fi(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedFi0(1:3,1:3,1:myNgrains,i,e) ! ...intermediate def grads
|
|
|
|
crystallite_Li(1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedLi0(1:3,1:3,1:myNgrains,i,e) ! ...intermediate velocity grads
|
|
|
|
crystallite_dPdF(1:3,1:3,1:3,1:3,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partioneddPdF0(1:3,1:3,1:3,1:3,1:myNgrains,i,e) ! ...stiffness
|
|
|
|
crystallite_Tstar_v(1:6,1:myNgrains,i,e) = &
|
|
|
|
crystallite_partionedTstar0_v(1:6,1:myNgrains,i,e) ! ...2nd PK stress
|
|
|
|
do g = 1, myNgrains
|
2016-01-15 05:49:44 +05:30
|
|
|
plasticState (phaseAt(g,i,e))%state( :,phasememberAt(g,i,e)) = &
|
|
|
|
plasticState (phaseAt(g,i,e))%partionedState0(:,phasememberAt(g,i,e))
|
|
|
|
do mySource = 1_pInt, phase_Nsources(phaseAt(g,i,e))
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%state( :,phasememberAt(g,i,e)) = &
|
|
|
|
sourceState(phaseAt(g,i,e))%p(mySource)%partionedState0(:,phasememberAt(g,i,e))
|
2015-05-28 22:32:23 +05:30
|
|
|
enddo
|
2015-10-15 00:06:19 +05:30
|
|
|
enddo
|
2015-05-28 22:32:23 +05:30
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
homogState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) ! ...internal homogenization state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
thermalState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) ! ...internal thermal state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
damageState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e)) ! ...internal damage state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
vacancyfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e))! ...internal vacancy transport state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
porosityState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e))! ...internal porosity state
|
|
|
|
forall(i = FEsolving_execIP(1,e):FEsolving_execIP(2,e), &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%sizeState > 0_pInt) &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%State( :,mappingHomogenization(1,i,e)) = &
|
|
|
|
hydrogenfluxState(mappingHomogenization(2,i,e))%subState0(:,mappingHomogenization(1,i,e))! ...internal hydrogen transport state
|
2015-10-15 00:06:19 +05:30
|
|
|
endif
|
2013-01-29 15:58:01 +05:30
|
|
|
endif converged
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2012-12-14 20:00:08 +05:30
|
|
|
if (materialpoint_subStep(i,e) > subStepMinHomog) then
|
|
|
|
materialpoint_requested(i,e) = .true.
|
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
|
|
|
materialpoint_subF(1:3,1:3,i,e) = materialpoint_subF0(1:3,1:3,i,e) + &
|
|
|
|
materialpoint_subStep(i,e) * (materialpoint_F(1:3,1:3,i,e) - materialpoint_F0(1:3,1:3,i,e))
|
|
|
|
materialpoint_subdt(i,e) = materialpoint_subStep(i,e) * dt
|
2013-01-29 15:58:01 +05:30
|
|
|
materialpoint_doneAndHappy(1:2,i,e) = [.false.,.true.]
|
2009-05-07 21:57:36 +05:30
|
|
|
endif
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo IpLooping1
|
|
|
|
enddo elementLooping1
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END PARALLEL DO
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2009-08-11 22:01:57 +05:30
|
|
|
NiterationMPstate = 0_pInt
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
convergenceLooping: do while (.not. terminallyIll .and. &
|
2010-09-02 02:34:02 +05:30
|
|
|
any( materialpoint_requested(:,FEsolving_execELem(1):FEsolving_execElem(2)) &
|
2009-05-07 21:57:36 +05:30
|
|
|
.and. .not. materialpoint_doneAndHappy(1,:,FEsolving_execELem(1):FEsolving_execElem(2)) &
|
2010-09-02 02:34:02 +05:30
|
|
|
) .and. &
|
2013-01-29 15:58:01 +05:30
|
|
|
NiterationMPstate < nMPstate)
|
2009-08-11 22:01:57 +05:30
|
|
|
NiterationMPstate = NiterationMPstate + 1
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! deformation partitioning
|
2015-10-15 00:06:19 +05:30
|
|
|
! based on materialpoint_subF0,.._subF,crystallite_partionedF0, and homogenization_state,
|
2009-05-07 21:57:36 +05:30
|
|
|
! results in crystallite_partionedF
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP PARALLEL DO PRIVATE(myNgrains)
|
2013-01-29 15:58:01 +05:30
|
|
|
elementLooping2: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
2009-05-07 21:57:36 +05:30
|
|
|
myNgrains = homogenization_Ngrains(mesh_element(3,e))
|
2013-01-29 15:58:01 +05:30
|
|
|
IpLooping2: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
|
|
|
if ( materialpoint_requested(i,e) .and. & ! process requested but...
|
|
|
|
.not. materialpoint_doneAndHappy(1,i,e)) then ! ...not yet done material points
|
|
|
|
call homogenization_partitionDeformation(i,e) ! partition deformation onto constituents
|
|
|
|
crystallite_dt(1:myNgrains,i,e) = materialpoint_subdt(i,e) ! propagate materialpoint dt to grains
|
|
|
|
crystallite_requested(1:myNgrains,i,e) = .true. ! request calculation for constituents
|
2009-08-27 17:40:06 +05:30
|
|
|
else
|
2013-01-29 15:58:01 +05:30
|
|
|
crystallite_requested(1:myNgrains,i,e) = .false. ! calculation for constituents not required anymore
|
2009-05-07 21:57:36 +05:30
|
|
|
endif
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo IpLooping2
|
|
|
|
enddo elementLooping2
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END PARALLEL DO
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! crystallite integration
|
2009-05-07 21:57:36 +05:30
|
|
|
! based on crystallite_partionedF0,.._partionedF
|
|
|
|
! incrementing by crystallite_dt
|
2014-08-08 02:38:34 +05:30
|
|
|
call crystallite_stressAndItsTangent(updateJaco) ! request stress and tangent calculation for constituent grains
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! state update
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP PARALLEL DO
|
2013-01-29 15:58:01 +05:30
|
|
|
elementLooping3: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
|
|
|
IpLooping3: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
2009-05-07 21:57:36 +05:30
|
|
|
if ( materialpoint_requested(i,e) .and. &
|
|
|
|
.not. materialpoint_doneAndHappy(1,i,e)) then
|
2009-08-11 22:01:57 +05:30
|
|
|
if (.not. all(crystallite_converged(:,i,e))) then
|
2013-01-29 15:58:01 +05:30
|
|
|
materialpoint_doneAndHappy(1:2,i,e) = [.true.,.false.]
|
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
|
|
|
materialpoint_converged(i,e) = .false.
|
2009-08-11 22:01:57 +05:30
|
|
|
else
|
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
|
|
|
materialpoint_doneAndHappy(1:2,i,e) = homogenization_updateState(i,e)
|
2015-05-28 22:32:23 +05:30
|
|
|
materialpoint_converged(i,e) = all(materialpoint_doneAndHappy(1:2,i,e)) ! converged if done and happy
|
2009-08-11 22:01:57 +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
|
|
|
!$OMP FLUSH(materialpoint_converged)
|
2010-11-03 22:52:48 +05:30
|
|
|
if (materialpoint_converged(i,e)) then
|
2012-07-05 15:24:50 +05:30
|
|
|
if (iand(debug_level(debug_homogenization), debug_levelBasic) /= 0_pInt) then
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP CRITICAL (distributionMPState)
|
|
|
|
debug_MaterialpointStateLoopdistribution(NiterationMPstate) = &
|
2013-01-29 15:58:01 +05:30
|
|
|
debug_MaterialpointStateLoopdistribution(NiterationMPstate) + 1_pInt
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP END CRITICAL (distributionMPState)
|
|
|
|
endif
|
2010-11-03 22:52:48 +05:30
|
|
|
endif
|
2009-05-07 21:57:36 +05:30
|
|
|
endif
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo IpLooping3
|
|
|
|
enddo elementLooping3
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END PARALLEL DO
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo convergenceLooping
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2010-10-01 17:48:49 +05:30
|
|
|
NiterationHomog = NiterationHomog + 1_pInt
|
2009-08-11 22:01:57 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo cutBackLooping
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2015-10-15 00:06:19 +05:30
|
|
|
if (.not. terminallyIll ) then
|
2013-01-29 15:58:01 +05:30
|
|
|
call crystallite_orientations() ! calculate crystal orientations
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP PARALLEL DO
|
2013-01-29 15:58:01 +05:30
|
|
|
elementLooping4: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
|
|
|
IpLooping4: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
2009-08-11 22:01:57 +05:30
|
|
|
call homogenization_averageStressAndItsTangent(i,e)
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo IpLooping4
|
|
|
|
enddo elementLooping4
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END PARALLEL DO
|
2010-09-02 02:34:02 +05:30
|
|
|
else
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP CRITICAL (write2out)
|
2013-01-29 15:58:01 +05:30
|
|
|
write(6,'(/,a,/)') '<< HOMOG >> Material Point terminally ill'
|
2011-03-21 16:01:17 +05:30
|
|
|
!$OMP END CRITICAL (write2out)
|
2010-03-19 19:44:08 +05:30
|
|
|
endif
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end subroutine materialpoint_stressAndItsTangent
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief parallelized calculation of result array at material points
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-19 00:27:28 +05:30
|
|
|
subroutine materialpoint_postResults
|
2013-01-29 15:58:01 +05:30
|
|
|
use FEsolving, only: &
|
|
|
|
FEsolving_execElem, &
|
|
|
|
FEsolving_execIP
|
|
|
|
use mesh, only: &
|
|
|
|
mesh_element
|
|
|
|
use material, only: &
|
2014-09-26 16:02:56 +05:30
|
|
|
mappingHomogenization, &
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
2016-01-15 05:49:44 +05:30
|
|
|
phaseAt, phasememberAt, &
|
2015-03-30 15:15:10 +05:30
|
|
|
homogenization_maxNgrains, &
|
|
|
|
material_Ncrystallite, &
|
2015-10-15 00:06:19 +05:30
|
|
|
material_Nphase, &
|
2015-05-28 22:32:23 +05:30
|
|
|
#else
|
2014-09-26 16:02:56 +05:30
|
|
|
homogState, &
|
2014-06-30 20:17:30 +05:30
|
|
|
thermalState, &
|
2015-05-28 22:32:23 +05:30
|
|
|
damageState, &
|
|
|
|
vacancyfluxState, &
|
|
|
|
porosityState, &
|
|
|
|
hydrogenfluxState, &
|
|
|
|
#endif
|
|
|
|
plasticState, &
|
|
|
|
sourceState, &
|
2014-06-30 20:17:30 +05:30
|
|
|
material_phase, &
|
2013-01-29 15:58:01 +05:30
|
|
|
homogenization_Ngrains, &
|
|
|
|
microstructure_crystallite
|
|
|
|
use constitutive, only: &
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
2015-05-28 22:32:23 +05:30
|
|
|
constitutive_plasticity_maxSizePostResults, &
|
|
|
|
constitutive_source_maxSizePostResults, &
|
2015-03-30 15:15:10 +05:30
|
|
|
#endif
|
2013-01-29 15:58:01 +05:30
|
|
|
constitutive_postResults
|
|
|
|
use crystallite, only: &
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
|
|
|
crystallite_maxSizePostResults, &
|
|
|
|
#endif
|
2013-01-29 15:58:01 +05:30
|
|
|
crystallite_sizePostResults, &
|
|
|
|
crystallite_postResults
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
implicit none
|
2013-04-26 18:53:36 +05:30
|
|
|
integer(pInt) :: &
|
2013-01-29 15:58:01 +05:30
|
|
|
thePos, &
|
|
|
|
theSize, &
|
|
|
|
myNgrains, &
|
|
|
|
myCrystallite, &
|
|
|
|
g, & !< grain number
|
|
|
|
i, & !< integration point number
|
|
|
|
e !< element number
|
2015-03-30 15:15:10 +05:30
|
|
|
#ifdef FEM
|
|
|
|
integer(pInt) :: &
|
|
|
|
myHomog, &
|
|
|
|
myPhase, &
|
|
|
|
crystalliteCtr(material_Ncrystallite, homogenization_maxNgrains), &
|
2015-10-15 00:06:19 +05:30
|
|
|
phaseCtr (material_Nphase, homogenization_maxNgrains)
|
2015-03-30 15:15:10 +05:30
|
|
|
real(pReal), dimension(1+crystallite_maxSizePostResults + &
|
2015-05-28 22:32:23 +05:30
|
|
|
1+constitutive_plasticity_maxSizePostResults + &
|
2015-10-15 00:06:19 +05:30
|
|
|
constitutive_source_maxSizePostResults) :: &
|
2015-03-30 15:15:10 +05:30
|
|
|
crystalliteResults
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
crystalliteCtr = 0_pInt; phaseCtr = 0_pInt
|
|
|
|
elementLooping: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
|
|
|
myNgrains = homogenization_Ngrains(mesh_element(3,e))
|
|
|
|
myCrystallite = microstructure_crystallite(mesh_element(4,e))
|
|
|
|
IpLooping: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
|
|
|
myHomog = mappingHomogenization(2,i,e)
|
|
|
|
thePos = mappingHomogenization(1,i,e)
|
|
|
|
homogOutput(myHomog)%output(1: &
|
2015-05-28 22:32:23 +05:30
|
|
|
homogOutput(myHomog)%sizeResults, &
|
2015-03-30 15:15:10 +05:30
|
|
|
thePos) = homogenization_postResults(i,e)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-03-30 15:15:10 +05:30
|
|
|
grainLooping :do g = 1,myNgrains
|
2016-01-15 05:49:44 +05:30
|
|
|
myPhase = phaseAt(g,i,e)
|
2015-03-30 15:15:10 +05:30
|
|
|
crystalliteResults(1:1+crystallite_sizePostResults(myCrystallite) + &
|
|
|
|
1+plasticState(myPhase)%sizePostResults + &
|
2015-05-28 22:32:23 +05:30
|
|
|
sum(sourceState(myPhase)%p(:)%sizePostResults)) = crystallite_postResults(g,i,e)
|
2015-03-30 15:15:10 +05:30
|
|
|
if (microstructure_crystallite(mesh_element(4,e)) == myCrystallite .and. &
|
|
|
|
homogenization_Ngrains (mesh_element(3,e)) >= g) then
|
|
|
|
crystalliteCtr(myCrystallite,g) = crystalliteCtr(myCrystallite,g) + 1_pInt
|
|
|
|
crystalliteOutput(myCrystallite,g)% &
|
2015-10-15 00:06:19 +05:30
|
|
|
output(1:crystalliteOutput(myCrystallite,g)%sizeResults,crystalliteCtr(myCrystallite,g)) = &
|
2015-05-28 22:32:23 +05:30
|
|
|
crystalliteResults(2:1+crystalliteOutput(myCrystallite,g)%sizeResults)
|
2015-03-30 15:15:10 +05:30
|
|
|
endif
|
|
|
|
if (material_phase(g,i,e) == myPhase) then
|
|
|
|
phaseCtr(myPhase,g) = phaseCtr(myPhase,g) + 1_pInt
|
|
|
|
phaseOutput(myPhase,g)% &
|
2015-10-15 00:06:19 +05:30
|
|
|
output(1:phaseOutput(myPhase,g)%sizeResults,phaseCtr(myPhase,g)) = &
|
2015-05-28 22:32:23 +05:30
|
|
|
crystalliteResults(3 + crystalliteOutput(myCrystallite,g)%sizeResults: &
|
2015-10-15 00:06:19 +05:30
|
|
|
1 + crystalliteOutput(myCrystallite,g)%sizeResults + &
|
|
|
|
1 + plasticState (myphase)%sizePostResults + &
|
2015-05-28 22:32:23 +05:30
|
|
|
sum(sourceState(myphase)%p(:)%sizePostResults))
|
2015-03-30 15:15:10 +05:30
|
|
|
endif
|
|
|
|
enddo grainLooping
|
|
|
|
enddo IpLooping
|
|
|
|
enddo elementLooping
|
|
|
|
#else
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2011-08-01 23:40:55 +05:30
|
|
|
!$OMP PARALLEL DO PRIVATE(myNgrains,myCrystallite,thePos,theSize)
|
2013-01-29 15:58:01 +05:30
|
|
|
elementLooping: do e = FEsolving_execElem(1),FEsolving_execElem(2)
|
2009-05-07 21:57:36 +05:30
|
|
|
myNgrains = homogenization_Ngrains(mesh_element(3,e))
|
2010-02-25 23:09:11 +05:30
|
|
|
myCrystallite = microstructure_crystallite(mesh_element(4,e))
|
2013-01-29 15:58:01 +05:30
|
|
|
IpLooping: do i = FEsolving_execIP(1,e),FEsolving_execIP(2,e)
|
2011-08-01 23:40:55 +05:30
|
|
|
thePos = 0_pInt
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
theSize = homogState (mappingHomogenization(2,i,e))%sizePostResults &
|
|
|
|
+ thermalState (mappingHomogenization(2,i,e))%sizePostResults &
|
|
|
|
+ damageState (mappingHomogenization(2,i,e))%sizePostResults &
|
|
|
|
+ vacancyfluxState (mappingHomogenization(2,i,e))%sizePostResults &
|
|
|
|
+ porosityState (mappingHomogenization(2,i,e))%sizePostResults &
|
|
|
|
+ hydrogenfluxState(mappingHomogenization(2,i,e))%sizePostResults
|
2012-08-10 21:28:17 +05:30
|
|
|
materialpoint_results(thePos+1,i,e) = real(theSize,pReal) ! tell size of homogenization results
|
2011-08-01 23:40:55 +05:30
|
|
|
thePos = thePos + 1_pInt
|
2011-03-29 12:57:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
if (theSize > 0_pInt) then ! any homogenization results to mention?
|
|
|
|
materialpoint_results(thePos+1:thePos+theSize,i,e) = homogenization_postResults(i,e) ! tell homogenization results
|
2011-08-01 23:40:55 +05:30
|
|
|
thePos = thePos + theSize
|
2009-05-07 21:57:36 +05:30
|
|
|
endif
|
2014-09-26 16:04:36 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
materialpoint_results(thePos+1,i,e) = real(myNgrains,pReal) ! tell number of grains at materialpoint
|
2011-11-23 14:39:00 +05:30
|
|
|
thePos = thePos + 1_pInt
|
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
grainLooping :do g = 1,myNgrains
|
2014-09-26 20:46:10 +05:30
|
|
|
theSize = 1 + crystallite_sizePostResults(myCrystallite) + &
|
2015-05-28 22:32:23 +05:30
|
|
|
1 + plasticState (material_phase(g,i,e))%sizePostResults + & !ToDo
|
|
|
|
sum(sourceState(material_phase(g,i,e))%p(:)%sizePostResults)
|
2014-03-12 13:03:51 +05:30
|
|
|
materialpoint_results(thePos+1:thePos+theSize,i,e) = crystallite_postResults(g,i,e) ! tell crystallite results
|
2011-08-01 23:40:55 +05:30
|
|
|
thePos = thePos + theSize
|
2013-01-29 15:58:01 +05:30
|
|
|
enddo grainLooping
|
|
|
|
enddo IpLooping
|
|
|
|
enddo elementLooping
|
2010-11-03 22:52:48 +05:30
|
|
|
!$OMP END PARALLEL DO
|
2015-03-30 15:15:10 +05:30
|
|
|
#endif
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end subroutine materialpoint_postResults
|
2015-10-15 00:06:19 +05:30
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief partition material point def grad onto constituents
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-11 21:31:53 +05:30
|
|
|
subroutine homogenization_partitionDeformation(ip,el)
|
2013-01-29 15:58:01 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_element
|
|
|
|
use material, only: &
|
|
|
|
homogenization_type, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
2014-03-14 04:50:50 +05:30
|
|
|
HOMOGENIZATION_NONE_ID, &
|
2013-11-27 13:34:05 +05:30
|
|
|
HOMOGENIZATION_ISOSTRAIN_ID, &
|
|
|
|
HOMOGENIZATION_RGC_ID
|
2013-01-29 15:58:01 +05:30
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_partionedF
|
|
|
|
use homogenization_isostrain, only: &
|
|
|
|
homogenization_isostrain_partitionDeformation
|
|
|
|
use homogenization_RGC, only: &
|
|
|
|
homogenization_RGC_partitionDeformation
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
implicit none
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-11 21:31:53 +05:30
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
2014-03-14 04:50:50 +05:30
|
|
|
|
2013-10-11 21:31:53 +05:30
|
|
|
chosenHomogenization: select case(homogenization_type(mesh_element(3,el)))
|
2014-03-14 04:50:50 +05:30
|
|
|
|
|
|
|
case (HOMOGENIZATION_NONE_ID) chosenHomogenization
|
|
|
|
crystallite_partionedF(1:3,1:3,1:homogenization_maxNgrains,ip,el) = 0.0_pReal
|
|
|
|
crystallite_partionedF(1:3,1:3,1:1,ip,el) = &
|
|
|
|
spread(materialpoint_subF(1:3,1:3,ip,el),3,1)
|
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_ISOSTRAIN_ID) chosenHomogenization
|
2012-08-10 21:28:17 +05:30
|
|
|
call homogenization_isostrain_partitionDeformation(&
|
2013-10-11 21:31:53 +05:30
|
|
|
crystallite_partionedF(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
materialpoint_subF(1:3,1:3,ip,el),&
|
|
|
|
el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_RGC_ID) chosenHomogenization
|
2014-08-21 23:18:20 +05:30
|
|
|
call homogenization_RGC_partitionDeformation(&
|
|
|
|
crystallite_partionedF(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
materialpoint_subF(1:3,1:3,ip,el),&
|
|
|
|
ip, &
|
|
|
|
el)
|
2013-01-29 15:58:01 +05:30
|
|
|
end select chosenHomogenization
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end subroutine homogenization_partitionDeformation
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-10-15 00:06:19 +05:30
|
|
|
!> @brief update the internal state of the homogenization scheme and tell whether "done" and
|
2012-08-10 21:28:17 +05:30
|
|
|
!> "happy" with result
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-11 21:31:53 +05:30
|
|
|
function homogenization_updateState(ip,el)
|
2013-01-29 15:58:01 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_element
|
|
|
|
use material, only: &
|
|
|
|
homogenization_type, &
|
2015-05-28 22:32:23 +05:30
|
|
|
thermal_type, &
|
|
|
|
damage_type, &
|
|
|
|
vacancyflux_type, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
2015-05-28 22:32:23 +05:30
|
|
|
HOMOGENIZATION_RGC_ID, &
|
|
|
|
THERMAL_adiabatic_ID, &
|
|
|
|
DAMAGE_local_ID, &
|
|
|
|
VACANCYFLUX_isochempot_ID
|
2013-01-29 15:58:01 +05:30
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_P, &
|
|
|
|
crystallite_dPdF, &
|
|
|
|
crystallite_partionedF,&
|
|
|
|
crystallite_partionedF0
|
|
|
|
use homogenization_RGC, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_RGC_updateState
|
2015-05-28 22:32:23 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_updateState
|
|
|
|
use damage_local, only: &
|
|
|
|
damage_local_updateState
|
|
|
|
use vacancyflux_isochempot, only: &
|
|
|
|
vacancyflux_isochempot_updateState
|
2013-11-27 13:34:05 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
implicit none
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-11 21:31:53 +05:30
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
2009-05-07 21:57:36 +05:30
|
|
|
logical, dimension(2) :: homogenization_updateState
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
homogenization_updateState = .true.
|
2013-10-11 21:31:53 +05:30
|
|
|
chosenHomogenization: select case(homogenization_type(mesh_element(3,el)))
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_RGC_ID) chosenHomogenization
|
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
|
|
|
homogenization_updateState = &
|
2015-05-28 22:32:23 +05:30
|
|
|
homogenization_updateState .and. &
|
2014-08-21 23:18:20 +05:30
|
|
|
homogenization_RGC_updateState(crystallite_P(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
2015-05-28 22:32:23 +05:30
|
|
|
crystallite_partionedF(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
crystallite_partionedF0(1:3,1:3,1:homogenization_maxNgrains,ip,el),&
|
|
|
|
materialpoint_subF(1:3,1:3,ip,el),&
|
|
|
|
materialpoint_subdt(ip,el), &
|
|
|
|
crystallite_dPdF(1:3,1:3,1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
ip, &
|
|
|
|
el)
|
2013-01-29 15:58:01 +05:30
|
|
|
end select chosenHomogenization
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
chosenThermal: select case (thermal_type(mesh_element(3,el)))
|
|
|
|
case (THERMAL_adiabatic_ID) chosenThermal
|
|
|
|
homogenization_updateState = &
|
|
|
|
homogenization_updateState .and. &
|
|
|
|
thermal_adiabatic_updateState(materialpoint_subdt(ip,el), &
|
|
|
|
ip, &
|
|
|
|
el)
|
|
|
|
end select chosenThermal
|
|
|
|
|
|
|
|
chosenDamage: select case (damage_type(mesh_element(3,el)))
|
|
|
|
case (DAMAGE_local_ID) chosenDamage
|
|
|
|
homogenization_updateState = &
|
|
|
|
homogenization_updateState .and. &
|
|
|
|
damage_local_updateState(materialpoint_subdt(ip,el), &
|
|
|
|
ip, &
|
|
|
|
el)
|
|
|
|
end select chosenDamage
|
|
|
|
|
|
|
|
chosenVacancyflux: select case (vacancyflux_type(mesh_element(3,el)))
|
|
|
|
case (VACANCYFLUX_isochempot_ID) chosenVacancyflux
|
|
|
|
homogenization_updateState = &
|
|
|
|
homogenization_updateState .and. &
|
|
|
|
vacancyflux_isochempot_updateState(materialpoint_subdt(ip,el), &
|
|
|
|
ip, &
|
|
|
|
el)
|
|
|
|
end select chosenVacancyflux
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end function homogenization_updateState
|
2009-05-07 21:57:36 +05:30
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief derive average stress and stiffness from constituent quantities
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-11 21:31:53 +05:30
|
|
|
subroutine homogenization_averageStressAndItsTangent(ip,el)
|
2013-01-29 15:58:01 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_element
|
|
|
|
use material, only: &
|
|
|
|
homogenization_type, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_maxNgrains, &
|
2014-03-14 04:50:50 +05:30
|
|
|
HOMOGENIZATION_NONE_ID, &
|
2013-11-27 13:34:05 +05:30
|
|
|
HOMOGENIZATION_ISOSTRAIN_ID, &
|
|
|
|
HOMOGENIZATION_RGC_ID
|
2013-01-29 15:58:01 +05:30
|
|
|
use crystallite, only: &
|
|
|
|
crystallite_P,crystallite_dPdF
|
|
|
|
use homogenization_isostrain, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_isostrain_averageStressAndItsTangent
|
2013-01-29 15:58:01 +05:30
|
|
|
use homogenization_RGC, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_RGC_averageStressAndItsTangent
|
2009-07-22 21:37:19 +05:30
|
|
|
|
|
|
|
implicit none
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-11 21:31:53 +05:30
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-10-11 21:31:53 +05:30
|
|
|
chosenHomogenization: select case(homogenization_type(mesh_element(3,el)))
|
2014-03-14 04:50:50 +05:30
|
|
|
case (HOMOGENIZATION_NONE_ID) chosenHomogenization
|
|
|
|
materialpoint_P(1:3,1:3,ip,el) = sum(crystallite_P(1:3,1:3,1:1,ip,el),3)
|
|
|
|
materialpoint_dPdF(1:3,1:3,1:3,1:3,ip,el) &
|
|
|
|
= sum(crystallite_dPdF(1:3,1:3,1:3,1:3,1:1,ip,el),5)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_ISOSTRAIN_ID) chosenHomogenization
|
2013-10-16 18:34:59 +05:30
|
|
|
call homogenization_isostrain_averageStressAndItsTangent(&
|
|
|
|
materialpoint_P(1:3,1:3,ip,el), &
|
|
|
|
materialpoint_dPdF(1:3,1:3,1:3,1:3,ip,el),&
|
|
|
|
crystallite_P(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
crystallite_dPdF(1:3,1:3,1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
el)
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_RGC_ID) chosenHomogenization
|
2013-10-16 18:34:59 +05:30
|
|
|
call homogenization_RGC_averageStressAndItsTangent(&
|
|
|
|
materialpoint_P(1:3,1:3,ip,el), &
|
|
|
|
materialpoint_dPdF(1:3,1:3,1:3,1:3,ip,el),&
|
|
|
|
crystallite_P(1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
crystallite_dPdF(1:3,1:3,1:3,1:3,1:homogenization_maxNgrains,ip,el), &
|
|
|
|
el)
|
2013-01-29 15:58:01 +05:30
|
|
|
end select chosenHomogenization
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end subroutine homogenization_averageStressAndItsTangent
|
2009-07-22 21:37:19 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-10-15 00:06:19 +05:30
|
|
|
!> @brief return array of homogenization results for post file inclusion. call only,
|
2013-01-29 15:58:01 +05:30
|
|
|
!> if homogenization_sizePostResults(i,e) > 0 !!
|
2012-08-10 21:28:17 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-10-11 21:31:53 +05:30
|
|
|
function homogenization_postResults(ip,el)
|
2013-01-29 15:58:01 +05:30
|
|
|
use mesh, only: &
|
|
|
|
mesh_element
|
|
|
|
use material, only: &
|
2014-09-26 16:02:56 +05:30
|
|
|
mappingHomogenization, &
|
|
|
|
homogState, &
|
2015-05-28 22:32:23 +05:30
|
|
|
thermalState, &
|
|
|
|
damageState, &
|
|
|
|
vacancyfluxState, &
|
|
|
|
porosityState, &
|
|
|
|
hydrogenfluxState, &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_type, &
|
2015-05-28 22:32:23 +05:30
|
|
|
thermal_type, &
|
|
|
|
damage_type, &
|
|
|
|
vacancyflux_type, &
|
|
|
|
porosity_type, &
|
|
|
|
hydrogenflux_type, &
|
2014-03-14 04:50:50 +05:30
|
|
|
HOMOGENIZATION_NONE_ID, &
|
2013-11-27 13:34:05 +05:30
|
|
|
HOMOGENIZATION_ISOSTRAIN_ID, &
|
2015-05-28 22:32:23 +05:30
|
|
|
HOMOGENIZATION_RGC_ID, &
|
|
|
|
THERMAL_isothermal_ID, &
|
|
|
|
THERMAL_adiabatic_ID, &
|
|
|
|
THERMAL_conduction_ID, &
|
|
|
|
DAMAGE_none_ID, &
|
|
|
|
DAMAGE_local_ID, &
|
|
|
|
DAMAGE_nonlocal_ID, &
|
|
|
|
VACANCYFLUX_isoconc_ID, &
|
|
|
|
VACANCYFLUX_isochempot_ID, &
|
|
|
|
VACANCYFLUX_cahnhilliard_ID, &
|
|
|
|
POROSITY_none_ID, &
|
|
|
|
POROSITY_phasefield_ID, &
|
|
|
|
HYDROGENFLUX_isoconc_ID, &
|
|
|
|
HYDROGENFLUX_cahnhilliard_ID
|
2013-01-29 15:58:01 +05:30
|
|
|
use homogenization_isostrain, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_isostrain_postResults
|
2013-01-29 15:58:01 +05:30
|
|
|
use homogenization_RGC, only: &
|
2013-11-27 13:34:05 +05:30
|
|
|
homogenization_RGC_postResults
|
2015-05-28 22:32:23 +05:30
|
|
|
use thermal_adiabatic, only: &
|
|
|
|
thermal_adiabatic_postResults
|
|
|
|
use thermal_conduction, only: &
|
|
|
|
thermal_conduction_postResults
|
|
|
|
use damage_local, only: &
|
|
|
|
damage_local_postResults
|
|
|
|
use damage_nonlocal, only: &
|
|
|
|
damage_nonlocal_postResults
|
|
|
|
use vacancyflux_isochempot, only: &
|
|
|
|
vacancyflux_isochempot_postResults
|
|
|
|
use vacancyflux_cahnhilliard, only: &
|
|
|
|
vacancyflux_cahnhilliard_postResults
|
|
|
|
use porosity_phasefield, only: &
|
|
|
|
porosity_phasefield_postResults
|
|
|
|
use hydrogenflux_cahnhilliard, only: &
|
|
|
|
hydrogenflux_cahnhilliard_postResults
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2009-05-07 21:57:36 +05:30
|
|
|
implicit none
|
2013-01-29 15:58:01 +05:30
|
|
|
integer(pInt), intent(in) :: &
|
2013-10-11 21:31:53 +05:30
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
2015-05-28 22:32:23 +05:30
|
|
|
real(pReal), dimension( homogState (mappingHomogenization(2,ip,el))%sizePostResults &
|
|
|
|
+ thermalState (mappingHomogenization(2,ip,el))%sizePostResults &
|
|
|
|
+ damageState (mappingHomogenization(2,ip,el))%sizePostResults &
|
|
|
|
+ vacancyfluxState (mappingHomogenization(2,ip,el))%sizePostResults &
|
|
|
|
+ porosityState (mappingHomogenization(2,ip,el))%sizePostResults &
|
|
|
|
+ hydrogenfluxState(mappingHomogenization(2,ip,el))%sizePostResults) :: &
|
2014-09-26 16:02:56 +05:30
|
|
|
homogenization_postResults
|
2015-05-28 22:32:23 +05:30
|
|
|
integer(pInt) :: &
|
|
|
|
startPos, endPos
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2009-05-07 21:57:36 +05:30
|
|
|
homogenization_postResults = 0.0_pReal
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
startPos = 1_pInt
|
|
|
|
endPos = homogState(mappingHomogenization(2,ip,el))%sizePostResults
|
2013-10-11 21:31:53 +05:30
|
|
|
chosenHomogenization: select case (homogenization_type(mesh_element(3,el)))
|
2014-03-14 04:50:50 +05:30
|
|
|
case (HOMOGENIZATION_NONE_ID) chosenHomogenization
|
|
|
|
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_ISOSTRAIN_ID) chosenHomogenization
|
2015-05-28 22:32:23 +05:30
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
homogenization_isostrain_postResults(&
|
2013-10-19 00:27:28 +05:30
|
|
|
ip, &
|
|
|
|
el, &
|
|
|
|
materialpoint_P(1:3,1:3,ip,el), &
|
|
|
|
materialpoint_F(1:3,1:3,ip,el))
|
2013-11-27 13:34:05 +05:30
|
|
|
case (HOMOGENIZATION_RGC_ID) chosenHomogenization
|
2015-05-28 22:32:23 +05:30
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
homogenization_RGC_postResults(&
|
2014-08-21 23:18:20 +05:30
|
|
|
ip, &
|
|
|
|
el, &
|
|
|
|
materialpoint_P(1:3,1:3,ip,el), &
|
|
|
|
materialpoint_F(1:3,1:3,ip,el))
|
2013-01-29 15:58:01 +05:30
|
|
|
end select chosenHomogenization
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + thermalState(mappingHomogenization(2,ip,el))%sizePostResults
|
|
|
|
chosenThermal: select case (thermal_type(mesh_element(3,el)))
|
|
|
|
case (THERMAL_isothermal_ID) chosenThermal
|
|
|
|
|
|
|
|
case (THERMAL_adiabatic_ID) chosenThermal
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
thermal_adiabatic_postResults(ip, el)
|
|
|
|
case (THERMAL_conduction_ID) chosenThermal
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
thermal_conduction_postResults(ip, el)
|
|
|
|
end select chosenThermal
|
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + damageState(mappingHomogenization(2,ip,el))%sizePostResults
|
|
|
|
chosenDamage: select case (damage_type(mesh_element(3,el)))
|
|
|
|
case (DAMAGE_none_ID) chosenDamage
|
|
|
|
|
|
|
|
case (DAMAGE_local_ID) chosenDamage
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
damage_local_postResults(ip, el)
|
|
|
|
|
|
|
|
case (DAMAGE_nonlocal_ID) chosenDamage
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
damage_nonlocal_postResults(ip, el)
|
|
|
|
end select chosenDamage
|
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + vacancyfluxState(mappingHomogenization(2,ip,el))%sizePostResults
|
|
|
|
chosenVacancyflux: select case (vacancyflux_type(mesh_element(3,el)))
|
|
|
|
case (VACANCYFLUX_isoconc_ID) chosenVacancyflux
|
|
|
|
|
|
|
|
case (VACANCYFLUX_isochempot_ID) chosenVacancyflux
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
vacancyflux_isochempot_postResults(ip, el)
|
|
|
|
case (VACANCYFLUX_cahnhilliard_ID) chosenVacancyflux
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
vacancyflux_cahnhilliard_postResults(ip, el)
|
|
|
|
end select chosenVacancyflux
|
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + porosityState(mappingHomogenization(2,ip,el))%sizePostResults
|
|
|
|
chosenPorosity: select case (porosity_type(mesh_element(3,el)))
|
|
|
|
case (POROSITY_none_ID) chosenPorosity
|
|
|
|
|
|
|
|
case (POROSITY_phasefield_ID) chosenPorosity
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
porosity_phasefield_postResults(ip, el)
|
|
|
|
end select chosenPorosity
|
|
|
|
|
|
|
|
startPos = endPos + 1_pInt
|
|
|
|
endPos = endPos + hydrogenfluxState(mappingHomogenization(2,ip,el))%sizePostResults
|
|
|
|
chosenHydrogenflux: select case (hydrogenflux_type(mesh_element(3,el)))
|
|
|
|
case (HYDROGENFLUX_isoconc_ID) chosenHydrogenflux
|
|
|
|
|
|
|
|
case (HYDROGENFLUX_cahnhilliard_ID) chosenHydrogenflux
|
|
|
|
homogenization_postResults(startPos:endPos) = &
|
|
|
|
hydrogenflux_cahnhilliard_postResults(ip, el)
|
|
|
|
end select chosenHydrogenflux
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end function homogenization_postResults
|
2014-09-26 16:04:36 +05:30
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end module homogenization
|