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
|
2019-06-15 20:12:16 +05:30
|
|
|
use prec
|
|
|
|
use IO
|
|
|
|
use config
|
|
|
|
use math
|
|
|
|
use material
|
|
|
|
use constitutive
|
|
|
|
use discretization
|
|
|
|
use thermal_isothermal
|
|
|
|
use thermal_conduction
|
|
|
|
use damage_none
|
|
|
|
use damage_nonlocal
|
2020-12-30 04:44:48 +05:30
|
|
|
use HDF5_utilities
|
2019-06-15 20:12:16 +05:30
|
|
|
use results
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2019-06-15 20:12:16 +05:30
|
|
|
implicit none
|
|
|
|
private
|
2019-06-11 13:18:07 +05:30
|
|
|
|
2020-03-29 23:34:51 +05:30
|
|
|
logical, public :: &
|
|
|
|
terminallyIll = .false. !< at least one material point is terminally ill
|
2020-04-14 13:08:32 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! General variables for the homogenization at a material point
|
2020-12-30 16:30:47 +05:30
|
|
|
real(pReal), dimension(:), allocatable, public :: &
|
|
|
|
homogenization_T
|
2020-12-16 17:18:45 +05:30
|
|
|
real(pReal), dimension(:,:,:), allocatable, public :: &
|
2020-12-16 15:51:24 +05:30
|
|
|
homogenization_F0, & !< def grad of IP at start of FE increment
|
|
|
|
homogenization_F !< def grad of IP to be reached at end of FE increment
|
2020-12-16 17:18:45 +05:30
|
|
|
real(pReal), dimension(:,:,:), allocatable, public :: & !, protected :: & Issue with ifort
|
2020-12-16 15:51:24 +05:30
|
|
|
homogenization_P !< first P--K stress of IP
|
2020-12-16 17:18:45 +05:30
|
|
|
real(pReal), dimension(:,:,:,:,:), allocatable, public :: & !, protected :: &
|
2020-12-16 15:51:24 +05:30
|
|
|
homogenization_dPdF !< tangent of first P--K stress at IP
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-03-17 05:09:32 +05:30
|
|
|
type :: tNumerics
|
|
|
|
integer :: &
|
2020-04-14 13:08:32 +05:30
|
|
|
nMPstate !< materialpoint state loop limit
|
2020-03-17 05:09:32 +05:30
|
|
|
real(pReal) :: &
|
|
|
|
subStepMinHomog, & !< minimum (relative) size of sub-step allowed during cutback in homogenization
|
|
|
|
subStepSizeHomog, & !< size of first substep when cutback in homogenization
|
|
|
|
stepIncreaseHomog !< increase of next substep size when previous substep converged in homogenization
|
|
|
|
end type tNumerics
|
|
|
|
|
|
|
|
type(tNumerics) :: num
|
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
interface
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
module subroutine mech_init(num_homog)
|
2020-06-16 22:17:19 +05:30
|
|
|
class(tNode), pointer, intent(in) :: &
|
2020-12-16 15:51:24 +05:30
|
|
|
num_homog !< pointer to mechanical homogenization numerics data
|
|
|
|
end subroutine mech_init
|
|
|
|
|
2020-12-30 16:30:47 +05:30
|
|
|
module subroutine thermal_init
|
|
|
|
end subroutine thermal_init
|
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
module subroutine mech_partition(subF,ip,el)
|
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
subF
|
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
|
|
|
end subroutine mech_partition
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2021-01-08 02:45:18 +05:30
|
|
|
module subroutine thermal_partition(T,ip,el)
|
|
|
|
real(pReal), intent(in) :: T
|
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
|
|
|
end subroutine thermal_partition
|
|
|
|
|
2020-12-29 05:09:23 +05:30
|
|
|
module subroutine mech_homogenize(dt,ip,el)
|
|
|
|
real(pReal), intent(in) :: dt
|
2020-12-16 15:51:24 +05:30
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
|
|
|
end subroutine mech_homogenize
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
module subroutine mech_results(group_base,h)
|
|
|
|
character(len=*), intent(in) :: group_base
|
|
|
|
integer, intent(in) :: h
|
|
|
|
end subroutine mech_results
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-28 14:25:54 +05:30
|
|
|
module function mech_updateState(subdt,subF,ip,el) result(doneAndHappy)
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
subdt !< current time step
|
|
|
|
real(pReal), intent(in), dimension(3,3) :: &
|
|
|
|
subF
|
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point
|
|
|
|
el !< element number
|
|
|
|
logical, dimension(2) :: doneAndHappy
|
|
|
|
end function mech_updateState
|
2019-06-15 20:12:16 +05:30
|
|
|
|
|
|
|
end interface
|
2013-01-29 15:58:01 +05:30
|
|
|
|
2019-06-15 20:12:16 +05:30
|
|
|
public :: &
|
|
|
|
homogenization_init, &
|
|
|
|
materialpoint_stressAndItsTangent, &
|
2020-12-23 11:28:54 +05:30
|
|
|
homogenization_forward, &
|
2020-12-30 04:44:48 +05:30
|
|
|
homogenization_results, &
|
|
|
|
homogenization_restartRead, &
|
|
|
|
homogenization_restartWrite
|
2012-08-10 21:28:17 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief module initialization
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-07-24 20:27:29 +05:30
|
|
|
subroutine homogenization_init
|
2012-05-08 20:27:06 +05:30
|
|
|
|
2020-06-16 22:17:19 +05:30
|
|
|
class (tNode) , pointer :: &
|
|
|
|
num_homog, &
|
2020-12-23 18:40:26 +05:30
|
|
|
num_homogGeneric
|
2020-07-02 01:50:22 +05:30
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
print'(/,a)', ' <<<+- homogenization init -+>>>'; flush(IO_STDOUT)
|
|
|
|
|
2020-09-13 14:09:17 +05:30
|
|
|
num_homog => config_numerics%get('homogenization',defaultVal=emptyDict)
|
2020-06-16 22:17:19 +05:30
|
|
|
num_homogGeneric => num_homog%get('generic',defaultVal=emptyDict)
|
|
|
|
|
2020-06-26 15:52:33 +05:30
|
|
|
num%nMPstate = num_homogGeneric%get_asInt ('nMPstate', defaultVal=10)
|
2020-06-16 22:17:19 +05:30
|
|
|
num%subStepMinHomog = num_homogGeneric%get_asFloat('subStepMin', defaultVal=1.0e-3_pReal)
|
|
|
|
num%subStepSizeHomog = num_homogGeneric%get_asFloat('subStepSize', defaultVal=0.25_pReal)
|
|
|
|
num%stepIncreaseHomog = num_homogGeneric%get_asFloat('stepIncrease', defaultVal=1.5_pReal)
|
|
|
|
|
2020-03-17 05:09:32 +05:30
|
|
|
if (num%nMPstate < 1) call IO_error(301,ext_msg='nMPstate')
|
|
|
|
if (num%subStepMinHomog <= 0.0_pReal) call IO_error(301,ext_msg='subStepMinHomog')
|
|
|
|
if (num%subStepSizeHomog <= 0.0_pReal) call IO_error(301,ext_msg='subStepSizeHomog')
|
|
|
|
if (num%stepIncreaseHomog <= 0.0_pReal) call IO_error(301,ext_msg='stepIncreaseHomog')
|
|
|
|
|
2020-12-16 15:51:24 +05:30
|
|
|
|
|
|
|
call mech_init(num_homog)
|
2021-01-08 02:45:18 +05:30
|
|
|
call thermal_init()
|
2020-12-16 15:51:24 +05:30
|
|
|
|
2021-01-08 02:45:18 +05:30
|
|
|
if (any(thermal_type == THERMAL_isothermal_ID)) call thermal_isothermal_init(homogenization_T)
|
|
|
|
if (any(thermal_type == THERMAL_conduction_ID)) call thermal_conduction_init(homogenization_T)
|
2020-12-16 15:51:24 +05:30
|
|
|
|
|
|
|
if (any(damage_type == DAMAGE_none_ID)) call damage_none_init
|
|
|
|
if (any(damage_type == DAMAGE_nonlocal_ID)) call damage_nonlocal_init
|
|
|
|
|
|
|
|
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-12-27 16:18:02 +05:30
|
|
|
subroutine materialpoint_stressAndItsTangent(dt,FEsolving_execIP,FEsolving_execElem)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2019-06-15 20:12:16 +05:30
|
|
|
real(pReal), intent(in) :: dt !< time increment
|
2020-12-27 16:18:02 +05:30
|
|
|
integer, dimension(2), intent(in) :: FEsolving_execElem, FEsolving_execIP
|
2019-06-15 20:12:16 +05:30
|
|
|
integer :: &
|
|
|
|
NiterationMPstate, &
|
2020-12-23 18:54:44 +05:30
|
|
|
ip, & !< integration point number
|
|
|
|
el, & !< element number
|
2021-01-17 14:00:42 +05:30
|
|
|
myNgrains, co, ce, ho, me, ph
|
2020-12-23 21:44:26 +05:30
|
|
|
real(pReal) :: &
|
2020-04-15 16:39:05 +05:30
|
|
|
subFrac, &
|
|
|
|
subStep
|
2020-12-23 21:44:26 +05:30
|
|
|
logical :: &
|
2020-04-15 16:39:05 +05:30
|
|
|
converged
|
2020-12-23 21:44:26 +05:30
|
|
|
logical, dimension(2) :: &
|
2020-04-15 16:39:05 +05:30
|
|
|
doneAndHappy
|
2020-07-03 20:15:11 +05:30
|
|
|
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2020-12-28 14:57:52 +05:30
|
|
|
!$OMP PARALLEL DO PRIVATE(ce,me,ho,myNgrains,NiterationMPstate,subFrac,converged,subStep,doneAndHappy)
|
2020-12-23 18:54:44 +05:30
|
|
|
do el = FEsolving_execElem(1),FEsolving_execElem(2)
|
2020-12-27 16:18:02 +05:30
|
|
|
ho = material_homogenizationAt(el)
|
|
|
|
myNgrains = homogenization_Nconstituents(ho)
|
2020-12-23 21:44:26 +05:30
|
|
|
do ip = FEsolving_execIP(1),FEsolving_execIP(2)
|
2020-12-28 14:57:52 +05:30
|
|
|
me = material_homogenizationMemberAt(ip,el)
|
2020-12-23 19:07:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! initialize restoration points
|
2020-12-23 18:54:44 +05:30
|
|
|
call constitutive_initializeRestorationPoints(ip,el)
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-23 21:44:26 +05:30
|
|
|
subFrac = 0.0_pReal
|
|
|
|
converged = .false. ! pretend failed step ...
|
|
|
|
subStep = 1.0_pReal/num%subStepSizeHomog ! ... larger then the requested calculation
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-29 17:31:33 +05:30
|
|
|
if (homogState(ho)%sizeState > 0) homogState(ho)%subState0(:,me) = homogState(ho)%State0(:,me)
|
2021-01-08 11:40:38 +05:30
|
|
|
if (damageState_h(ho)%sizeState > 0) damageState_h(ho)%subState0(:,me) = damageState_h(ho)%State0(:,me)
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-24 17:14:26 +05:30
|
|
|
cutBackLooping: do while (.not. terminallyIll .and. subStep > num%subStepMinHomog)
|
2015-10-15 00:06:19 +05:30
|
|
|
|
2020-12-23 21:44:26 +05:30
|
|
|
if (converged) then
|
|
|
|
subFrac = subFrac + subStep
|
|
|
|
subStep = min(1.0_pReal-subFrac,num%stepIncreaseHomog*subStep) ! introduce flexibility for step increase/acceleration
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-23 21:44:26 +05:30
|
|
|
steppingNeeded: if (subStep > num%subStepMinHomog) then
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-04-15 12:23:25 +05:30
|
|
|
! wind forward grain starting point
|
2020-12-23 18:54:44 +05:30
|
|
|
call constitutive_windForward(ip,el)
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-29 17:31:33 +05:30
|
|
|
if(homogState(ho)%sizeState > 0) homogState(ho)%subState0(:,me) = homogState(ho)%State(:,me)
|
2021-01-08 11:40:38 +05:30
|
|
|
if(damageState_h(ho)%sizeState > 0) damageState_h(ho)%subState0(:,me) = damageState_h(ho)%State(:,me)
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2019-06-15 20:12:16 +05:30
|
|
|
endif steppingNeeded
|
2020-12-28 01:50:54 +05:30
|
|
|
elseif ( (myNgrains == 1 .and. subStep <= 1.0 ) .or. & ! single grain already tried internal subStepping in crystallite
|
2020-12-23 21:44:26 +05:30
|
|
|
num%subStepSizeHomog * subStep <= num%subStepMinHomog ) then ! would require too small subStep
|
2015-05-28 22:32:23 +05:30
|
|
|
! cutback makes no sense
|
2020-12-28 01:50:54 +05:30
|
|
|
if (.not. terminallyIll) & ! so first signals terminally ill...
|
|
|
|
print*, ' Integration point ', ip,' at element ', el, ' terminally ill'
|
|
|
|
terminallyIll = .true. ! ...and kills all others
|
|
|
|
else ! cutback makes sense
|
|
|
|
subStep = num%subStepSizeHomog * subStep ! crystallite had severe trouble, so do a significant cutback
|
|
|
|
|
2020-12-28 02:15:31 +05:30
|
|
|
call constitutive_restore(ip,el,subStep < 1.0_pReal)
|
2020-12-28 01:50:54 +05:30
|
|
|
|
2020-12-29 17:31:33 +05:30
|
|
|
if(homogState(ho)%sizeState > 0) homogState(ho)%State(:,me) = homogState(ho)%subState0(:,me)
|
2021-01-08 11:40:38 +05:30
|
|
|
if(damageState_h(ho)%sizeState > 0) damageState_h(ho)%State(:,me) = damageState_h(ho)%subState0(:,me)
|
2020-04-15 16:39:05 +05:30
|
|
|
endif
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-27 18:03:14 +05:30
|
|
|
if (subStep > num%subStepMinHomog) doneAndHappy = [.false.,.true.]
|
2019-06-15 20:12:16 +05:30
|
|
|
|
2020-12-23 16:55:56 +05:30
|
|
|
NiterationMPstate = 0
|
2021-01-17 14:00:42 +05:30
|
|
|
convergenceLooping: do while (.not. (terminallyIll .or. doneAndHappy(1)) &
|
|
|
|
.and. NiterationMPstate < num%nMPstate)
|
2020-12-23 16:55:56 +05:30
|
|
|
NiterationMPstate = NiterationMPstate + 1
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2013-01-29 15:58:01 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! deformation partitioning
|
2020-12-23 16:07:00 +05:30
|
|
|
|
2020-12-27 18:03:14 +05:30
|
|
|
if (.not. doneAndHappy(1)) then
|
2020-12-23 18:54:44 +05:30
|
|
|
ce = (el-1)*discretization_nIPs + ip
|
2021-01-17 14:00:42 +05:30
|
|
|
call mech_partition( homogenization_F0(1:3,1:3,ce) &
|
|
|
|
+ (homogenization_F(1:3,1:3,ce)-homogenization_F0(1:3,1:3,ce))*(subStep+subFrac), &
|
|
|
|
ip,el)
|
2020-12-23 21:44:26 +05:30
|
|
|
converged = .true.
|
2020-12-23 16:55:56 +05:30
|
|
|
do co = 1, myNgrains
|
2020-12-23 21:44:26 +05:30
|
|
|
converged = converged .and. crystallite_stress(dt*subStep,co,ip,el)
|
2020-12-23 16:55:56 +05:30
|
|
|
enddo
|
2009-05-07 21:57:36 +05:30
|
|
|
|
2020-12-23 21:44:26 +05:30
|
|
|
if (.not. converged) then
|
|
|
|
doneAndHappy = [.true.,.false.]
|
2019-06-15 20:12:16 +05:30
|
|
|
else
|
2020-12-23 18:54:44 +05:30
|
|
|
ce = (el-1)*discretization_nIPs + ip
|
2020-12-28 14:25:54 +05:30
|
|
|
doneAndHappy = mech_updateState(dt*subStep, &
|
|
|
|
homogenization_F0(1:3,1:3,ce) &
|
|
|
|
+ (homogenization_F(1:3,1:3,ce)-homogenization_F0(1:3,1:3,ce)) &
|
2020-12-23 21:44:26 +05:30
|
|
|
*(subStep+subFrac), &
|
2020-12-28 14:25:54 +05:30
|
|
|
ip,el)
|
2020-12-23 21:44:26 +05:30
|
|
|
converged = all(doneAndHappy)
|
2019-06-15 20:12:16 +05:30
|
|
|
endif
|
|
|
|
endif
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-23 16:55:56 +05:30
|
|
|
enddo convergenceLooping
|
2020-12-24 17:14:26 +05:30
|
|
|
enddo cutBackLooping
|
2020-12-23 19:07:12 +05:30
|
|
|
enddo
|
2020-12-24 17:14:26 +05:30
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
2020-12-23 19:07:12 +05:30
|
|
|
|
2019-06-15 20:12:16 +05:30
|
|
|
if (.not. terminallyIll ) then
|
2021-01-17 14:40:30 +05:30
|
|
|
!$OMP PARALLEL DO PRIVATE(ho,ph,ce)
|
2021-01-17 14:00:42 +05:30
|
|
|
do el = FEsolving_execElem(1),FEsolving_execElem(2)
|
|
|
|
if (terminallyIll) continue
|
|
|
|
ho = material_homogenizationAt(el)
|
|
|
|
do ip = FEsolving_execIP(1),FEsolving_execIP(2)
|
2021-01-17 14:40:30 +05:30
|
|
|
ce = (el-1)*discretization_nIPs + ip
|
|
|
|
call thermal_partition(homogenization_T(ce),ip,el)
|
2021-01-17 14:00:42 +05:30
|
|
|
do co = 1, homogenization_Nconstituents(ho)
|
|
|
|
ph = material_phaseAt(co,el)
|
|
|
|
call constitutive_thermal_initializeRestorationPoints(ph,material_phaseMemberAt(co,ip,el))
|
|
|
|
if (.not. thermal_stress(dt,ph,material_phaseMemberAt(co,ip,el))) then
|
|
|
|
if (.not. terminallyIll) & ! so first signals terminally ill...
|
|
|
|
print*, ' Integration point ', ip,' at element ', el, ' terminally ill'
|
|
|
|
terminallyIll = .true. ! ...and kills all others
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
|
|
|
|
!$OMP PARALLEL DO PRIVATE(ho)
|
2020-12-23 18:54:44 +05:30
|
|
|
elementLooping3: do el = FEsolving_execElem(1),FEsolving_execElem(2)
|
2020-12-27 16:18:02 +05:30
|
|
|
ho = material_homogenizationAt(el)
|
2020-12-23 18:54:44 +05:30
|
|
|
IpLooping3: do ip = FEsolving_execIP(1),FEsolving_execIP(2)
|
2021-01-17 14:00:42 +05:30
|
|
|
do co = 1, homogenization_Nconstituents(ho)
|
2020-12-27 16:18:02 +05:30
|
|
|
call crystallite_orientations(co,ip,el)
|
|
|
|
enddo
|
2020-12-29 05:09:23 +05:30
|
|
|
call mech_homogenize(dt,ip,el)
|
2020-12-23 14:35:02 +05:30
|
|
|
enddo IpLooping3
|
|
|
|
enddo elementLooping3
|
2019-06-15 20:12:16 +05:30
|
|
|
!$OMP END PARALLEL DO
|
|
|
|
else
|
2020-09-18 02:27:56 +05:30
|
|
|
print'(/,a,/)', ' << HOMOG >> Material Point terminally ill'
|
2019-06-15 20:12:16 +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
|
|
|
|
|
|
|
|
2019-04-30 22:15:16 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief writes homogenization results to HDF5 output file
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine homogenization_results
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-28 13:57:48 +05:30
|
|
|
integer :: ho
|
2020-11-06 05:29:12 +05:30
|
|
|
character(len=:), allocatable :: group_base,group
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-28 13:48:20 +05:30
|
|
|
|
2020-12-22 13:15:01 +05:30
|
|
|
call results_closeGroup(results_addGroup('current/homogenization/'))
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-28 13:57:48 +05:30
|
|
|
do ho=1,size(material_name_homogenization)
|
|
|
|
group_base = 'current/homogenization/'//trim(material_name_homogenization(ho))
|
2019-12-10 11:44:39 +05:30
|
|
|
call results_closeGroup(results_addGroup(group_base))
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2020-12-28 13:57:48 +05:30
|
|
|
call mech_results(group_base,ho)
|
2019-07-16 05:38:18 +05:30
|
|
|
|
2019-12-10 11:44:39 +05:30
|
|
|
group = trim(group_base)//'/damage'
|
|
|
|
call results_closeGroup(results_addGroup(group))
|
2020-12-28 13:57:48 +05:30
|
|
|
select case(damage_type(ho))
|
2019-12-10 11:44:39 +05:30
|
|
|
case(DAMAGE_NONLOCAL_ID)
|
2020-12-28 13:57:48 +05:30
|
|
|
call damage_nonlocal_results(ho,group)
|
2019-12-10 11:44:39 +05:30
|
|
|
end select
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2019-12-10 11:44:39 +05:30
|
|
|
group = trim(group_base)//'/thermal'
|
|
|
|
call results_closeGroup(results_addGroup(group))
|
2020-12-28 13:57:48 +05:30
|
|
|
select case(thermal_type(ho))
|
2019-12-11 00:55:19 +05:30
|
|
|
case(THERMAL_CONDUCTION_ID)
|
2020-12-28 13:57:48 +05:30
|
|
|
call thermal_conduction_results(ho,group)
|
2019-12-11 00:55:19 +05:30
|
|
|
end select
|
2020-03-17 05:09:32 +05:30
|
|
|
|
2019-12-10 11:44:39 +05:30
|
|
|
enddo
|
2019-12-19 00:35:51 +05:30
|
|
|
|
2019-04-30 22:15:16 +05:30
|
|
|
end subroutine homogenization_results
|
|
|
|
|
2020-12-23 11:28:54 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief Forward data after successful increment.
|
|
|
|
! ToDo: Any guessing for the current states possible?
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine homogenization_forward
|
|
|
|
|
2020-12-23 14:35:02 +05:30
|
|
|
integer :: ho
|
2020-12-23 11:28:54 +05:30
|
|
|
|
2020-12-28 13:48:20 +05:30
|
|
|
|
2020-12-23 11:28:54 +05:30
|
|
|
do ho = 1, size(material_name_homogenization)
|
|
|
|
homogState (ho)%state0 = homogState (ho)%state
|
2021-01-08 11:40:38 +05:30
|
|
|
damageState_h(ho)%state0 = damageState_h(ho)%state
|
2020-12-23 11:28:54 +05:30
|
|
|
enddo
|
|
|
|
|
|
|
|
end subroutine homogenization_forward
|
|
|
|
|
2020-12-30 04:44:48 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine homogenization_restartWrite(fileHandle)
|
|
|
|
|
|
|
|
integer(HID_T), intent(in) :: fileHandle
|
|
|
|
|
|
|
|
integer(HID_T), dimension(2) :: groupHandle
|
|
|
|
integer :: ho
|
|
|
|
|
|
|
|
|
|
|
|
groupHandle(1) = HDF5_addGroup(fileHandle,'homogenization')
|
|
|
|
|
|
|
|
do ho = 1, size(material_name_homogenization)
|
|
|
|
|
|
|
|
groupHandle(2) = HDF5_addGroup(groupHandle(1),material_name_homogenization(ho))
|
|
|
|
|
|
|
|
call HDF5_read(groupHandle(2),homogState(ho)%state,'omega') ! ToDo: should be done by mech
|
|
|
|
|
|
|
|
call HDF5_closeGroup(groupHandle(2))
|
|
|
|
|
|
|
|
enddo
|
|
|
|
|
|
|
|
call HDF5_closeGroup(groupHandle(1))
|
|
|
|
|
|
|
|
end subroutine homogenization_restartWrite
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine homogenization_restartRead(fileHandle)
|
|
|
|
|
|
|
|
integer(HID_T), intent(in) :: fileHandle
|
|
|
|
|
|
|
|
integer(HID_T), dimension(2) :: groupHandle
|
|
|
|
integer :: ho
|
|
|
|
|
|
|
|
|
|
|
|
groupHandle(1) = HDF5_openGroup(fileHandle,'homogenization')
|
|
|
|
|
|
|
|
do ho = 1, size(material_name_homogenization)
|
|
|
|
|
|
|
|
groupHandle(2) = HDF5_openGroup(groupHandle(1),material_name_homogenization(ho))
|
|
|
|
|
|
|
|
call HDF5_write(groupHandle(2),homogState(ho)%state,'omega') ! ToDo: should be done by mech
|
|
|
|
|
|
|
|
call HDF5_closeGroup(groupHandle(2))
|
|
|
|
|
|
|
|
enddo
|
|
|
|
|
|
|
|
call HDF5_closeGroup(groupHandle(1))
|
|
|
|
|
|
|
|
end subroutine homogenization_restartRead
|
|
|
|
|
|
|
|
|
2012-08-10 21:28:17 +05:30
|
|
|
end module homogenization
|