DAMASK_EICMD/src/phase.f90

678 lines
24 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
2013-03-06 20:11:15 +05:30
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
2020-07-15 18:05:21 +05:30
!> @brief elasticity, plasticity, damage & thermal internal microstructure state
!--------------------------------------------------------------------------------------------------
module phase
use prec
2021-11-25 10:51:56 +05:30
use constants
2019-06-15 19:10:22 +05:30
use math
use rotations
2022-01-31 19:35:15 +05:30
use polynomials
2022-12-08 02:25:20 +05:30
use tables
2019-06-15 19:10:22 +05:30
use IO
use config
use material
2023-01-19 22:07:45 +05:30
use result
2023-07-15 23:58:57 +05:30
use crystal
2019-06-15 19:10:22 +05:30
use discretization
use parallelization
2021-07-08 17:57:04 +05:30
use HDF5
use HDF5_utilities
2021-01-27 11:40:53 +05:30
implicit none(type,external)
2019-06-15 19:10:22 +05:30
private
2020-12-23 22:00:19 +05:30
2022-03-06 12:58:48 +05:30
type :: tState
integer :: &
sizeState = 0, & !< size of state
sizeDotState = 0, & !< size of dot state, i.e. state(1:sizeDot) follows time evolution by dotState rates
offsetDeltaState = 0, & !< index offset of delta state
sizeDeltaState = 0 !< size of delta state, i.e. state(offset+1:offset+sizeDelta) follows time evolution by deltaState increments
real(pREAL), allocatable, dimension(:) :: &
2022-03-06 12:58:48 +05:30
atol
! http://stackoverflow.com/questions/3948210
real(pREAL), pointer, dimension(:,:), contiguous :: & !< is basically an allocatable+target, but in a type needs to be pointer
2022-03-06 12:58:48 +05:30
state0, &
state, & !< state
dotState, & !< rate of state change
deltaState !< increment of state change
real(pREAL), pointer, dimension(:,:) :: &
2022-03-06 12:58:48 +05:30
deltaState2
end type
type, extends(tState) :: tPlasticState
logical :: nonlocal = .false.
end type
type :: tSourceState
type(tState), dimension(:), allocatable :: p !< tState for each active source mechanism in a phase
end type
enum, bind(c); enumerator :: &
UNDEFINED, &
MECHANICAL_PLASTICITY_NONE, &
MECHANICAL_PLASTICITY_ISOTROPIC, &
MECHANICAL_PLASTICITY_PHENOPOWERLAW, &
MECHANICAL_PLASTICITY_KINEHARDENING, &
MECHANICAL_PLASTICITY_DISLOTWIN, &
MECHANICAL_PLASTICITY_DISLOTUNGSTEN, &
MECHANICAL_PLASTICITY_NONLOCAL, &
MECHANICAL_EIGEN_THERMALEXPANSION, &
DAMAGE_ISOBRITTLE, &
DAMAGE_ANISOBRITTLE, &
THERMAL_SOURCE_DISSIPATION, &
THERMAL_SOURCE_EXTERNALHEAT
end enum
integer(kind(UNDEFINED)), dimension(:), allocatable :: &
mechanical_plasticity_type, & !< plasticity of each phase
2023-07-24 22:00:02 +05:30
damage_type !< damage type of each phase
integer(kind(UNDEFINED)), dimension(:,:), allocatable :: &
2023-07-24 22:00:02 +05:30
thermal_source_type, &
mechanical_eigen_kinematics_type
character(len=2), allocatable, dimension(:) :: phase_lattice
real(pREAL), allocatable, dimension(:) :: phase_cOverA
real(pREAL), allocatable, dimension(:) :: phase_rho
2021-05-22 20:51:07 +05:30
type(tRotationContainer), dimension(:), allocatable :: &
2021-07-24 18:33:26 +05:30
phase_O_0, &
phase_O
2020-12-21 00:50:39 +05:30
type :: tNumerics
integer :: &
iJacoLpresiduum, & !< frequency of Jacobian update of residuum in Lp
iJacoLiresiduum, & !< frequency of Jacobian update of residuum in Li
nState, & !< state loop limit
2023-07-26 02:34:29 +05:30
nStress_Lp, & !< stress loop limit for Lp
nStress_Li !< stress loop limit for Li
real(pREAL) :: &
2023-07-25 00:33:35 +05:30
stepMinCryst, & !< minimum (relative) size of sub-step allowed during cutback
stepSizeCryst, & !< size of first substep when cutback
stepSizeLp, & !< size of first substep when cutback in Lp calculation
stepSizeLi, & !< size of first substep when cutback in Li calculation
stepIncreaseCryst, & !< increase of next substep size when previous substep converged
rtol_crystalliteState, &
rtol_Lp, & !< relative tolerance in stress loop for Lp
atol_Lp, & !< absolute tolerance in stress loop for Lp
rtol_Li, & !< relative tolerance in stress loop for Li
atol_Li !< absolute tolerance in stress loop for Li
end type tNumerics
type(tNumerics) :: num ! numerics parameters. Better name?
2020-08-15 19:32:10 +05:30
type(tPlasticState), allocatable, dimension(:), public :: &
plasticState
type(tState), allocatable, dimension(:), public :: &
2021-02-13 12:25:32 +05:30
damageState
2020-08-15 19:32:10 +05:30
interface
2020-12-21 16:44:09 +05:30
! == cleaned:begin =================================================================================
module subroutine mechanical_init(phases,num_mech)
type(tDict), pointer :: phases, num_mech
2021-02-09 03:51:53 +05:30
end subroutine mechanical_init
module subroutine damage_init
end subroutine damage_init
2020-07-03 20:15:11 +05:30
module subroutine thermal_init(phases)
type(tDict), pointer :: phases
end subroutine thermal_init
2020-12-21 16:44:09 +05:30
2023-01-19 22:07:45 +05:30
module subroutine mechanical_result(group,ph)
2020-12-21 16:44:09 +05:30
character(len=*), intent(in) :: group
integer, intent(in) :: ph
2023-01-19 22:07:45 +05:30
end subroutine mechanical_result
2020-12-21 16:44:09 +05:30
2023-01-19 22:07:45 +05:30
module subroutine damage_result(group,ph)
2020-12-22 16:50:00 +05:30
character(len=*), intent(in) :: group
integer, intent(in) :: ph
2023-01-19 22:07:45 +05:30
end subroutine damage_result
2020-12-22 16:50:00 +05:30
2023-01-19 22:07:45 +05:30
module subroutine thermal_result(group,ph)
2022-02-19 23:26:41 +05:30
character(len=*), intent(in) :: group
integer, intent(in) :: ph
2023-01-19 22:07:45 +05:30
end subroutine thermal_result
2022-02-19 23:26:41 +05:30
2021-02-09 03:51:53 +05:30
module subroutine mechanical_forward()
end subroutine mechanical_forward
2020-12-31 14:24:13 +05:30
2021-04-06 13:52:47 +05:30
module subroutine damage_forward()
end subroutine damage_forward
2020-12-31 14:24:13 +05:30
module subroutine thermal_forward()
end subroutine thermal_forward
2020-12-23 03:57:56 +05:30
2021-02-09 03:51:53 +05:30
module subroutine mechanical_restore(ce,includeL)
2021-01-24 14:09:40 +05:30
integer, intent(in) :: ce
2020-12-31 14:24:13 +05:30
logical, intent(in) :: includeL
2021-02-09 03:51:53 +05:30
end subroutine mechanical_restore
2020-12-28 02:15:31 +05:30
2021-07-17 03:02:08 +05:30
module subroutine damage_restore(ce)
integer, intent(in) :: ce
end subroutine damage_restore
2020-12-31 14:24:13 +05:30
2021-07-17 15:20:21 +05:30
module function phase_mechanical_dPdF(Delta_t,co,ce) result(dPdF)
real(pREAL), intent(in) :: Delta_t
2020-12-30 02:01:22 +05:30
integer, intent(in) :: &
co, & !< counter in constituent loop
ce
real(pREAL), dimension(3,3,3,3) :: dPdF
2021-02-09 03:51:53 +05:30
end function phase_mechanical_dPdF
2020-12-30 02:01:22 +05:30
2021-02-09 03:51:53 +05:30
module subroutine mechanical_restartWrite(groupHandle,ph)
2020-12-30 04:44:48 +05:30
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
2021-02-09 03:51:53 +05:30
end subroutine mechanical_restartWrite
2020-12-30 04:44:48 +05:30
2022-01-11 01:04:50 +05:30
module subroutine thermal_restartWrite(groupHandle,ph)
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
end subroutine thermal_restartWrite
module subroutine damage_restartWrite(groupHandle,ph)
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
end subroutine damage_restartWrite
2021-02-09 03:51:53 +05:30
module subroutine mechanical_restartRead(groupHandle,ph)
2020-12-30 04:44:48 +05:30
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
2021-02-09 03:51:53 +05:30
end subroutine mechanical_restartRead
2020-12-30 04:44:48 +05:30
2022-01-11 01:04:50 +05:30
module subroutine thermal_restartRead(groupHandle,ph)
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
end subroutine thermal_restartRead
2020-12-30 14:24:06 +05:30
module subroutine damage_restartRead(groupHandle,ph)
integer(HID_T), intent(in) :: groupHandle
integer, intent(in) :: ph
end subroutine damage_restartRead
2021-04-13 00:51:15 +05:30
module function mechanical_S(ph,en) result(S)
integer, intent(in) :: ph,en
real(pREAL), dimension(3,3) :: S
2021-02-09 03:51:53 +05:30
end function mechanical_S
2020-12-30 14:24:06 +05:30
2021-04-13 00:51:15 +05:30
module function mechanical_L_p(ph,en) result(L_p)
integer, intent(in) :: ph,en
real(pREAL), dimension(3,3) :: L_p
2021-02-09 03:51:53 +05:30
end function mechanical_L_p
2020-12-30 14:24:06 +05:30
2021-04-13 00:51:15 +05:30
module function mechanical_F_e(ph,en) result(F_e)
integer, intent(in) :: ph,en
real(pREAL), dimension(3,3) :: F_e
2021-02-09 03:51:53 +05:30
end function mechanical_F_e
2020-12-30 14:24:06 +05:30
module function mechanical_F_i(ph,en) result(F_i)
integer, intent(in) :: ph,en
real(pREAL), dimension(3,3) :: F_i
end function mechanical_F_i
2021-04-08 17:01:21 +05:30
module function phase_F(co,ce) result(F)
integer, intent(in) :: co, ce
real(pREAL), dimension(3,3) :: F
2021-04-08 17:01:21 +05:30
end function phase_F
2020-12-30 15:33:13 +05:30
2021-04-08 00:36:29 +05:30
module function phase_P(co,ce) result(P)
integer, intent(in) :: co, ce
real(pREAL), dimension(3,3) :: P
2021-04-08 00:36:29 +05:30
end function phase_P
2021-01-21 01:24:31 +05:30
pure module function thermal_T(ph,en) result(T)
2021-04-25 11:36:52 +05:30
integer, intent(in) :: ph,en
real(pREAL) :: T
end function thermal_T
2020-12-30 14:24:06 +05:30
2021-04-25 11:36:52 +05:30
module function thermal_dot_T(ph,en) result(dot_T)
integer, intent(in) :: ph,en
real(pREAL) :: dot_T
end function thermal_dot_T
module function damage_phi(ph,en) result(phi)
integer, intent(in) :: ph,en
real(pREAL) :: phi
2021-02-13 11:45:57 +05:30
end function damage_phi
2020-12-30 16:30:47 +05:30
2021-04-08 00:36:29 +05:30
module subroutine phase_set_F(F,co,ce)
real(pREAL), dimension(3,3), intent(in) :: F
2021-02-15 23:13:51 +05:30
integer, intent(in) :: co, ce
2021-04-08 00:36:29 +05:30
end subroutine phase_set_F
2020-12-30 14:24:06 +05:30
2021-02-09 03:51:53 +05:30
module subroutine phase_thermal_setField(T,dot_T, co,ce)
real(pREAL), intent(in) :: T, dot_T
2021-04-11 16:51:27 +05:30
integer, intent(in) :: co, ce
2021-02-09 03:51:53 +05:30
end subroutine phase_thermal_setField
2021-01-24 17:56:01 +05:30
2021-04-08 17:01:21 +05:30
module subroutine phase_set_phi(phi,co,ce)
real(pREAL), intent(in) :: phi
2021-01-21 01:24:31 +05:30
integer, intent(in) :: co, ce
2021-04-08 17:01:21 +05:30
end subroutine phase_set_phi
2020-12-30 16:30:47 +05:30
2021-04-11 16:51:27 +05:30
module function phase_mu_phi(co,ce) result(mu)
integer, intent(in) :: co, ce
real(pREAL) :: mu
2021-04-11 16:51:27 +05:30
end function phase_mu_phi
module function phase_K_phi(co,ce) result(K)
integer, intent(in) :: co, ce
real(pREAL), dimension(3,3) :: K
2021-04-11 16:51:27 +05:30
end function phase_K_phi
module function phase_mu_T(co,ce) result(mu)
integer, intent(in) :: co, ce
real(pREAL) :: mu
2021-04-11 16:51:27 +05:30
end function phase_mu_T
module function phase_K_T(co,ce) result(K)
integer, intent(in) :: co, ce
real(pREAL), dimension(3,3) :: K
2021-04-11 16:51:27 +05:30
end function phase_K_T
2020-12-30 16:30:47 +05:30
2020-12-23 11:28:54 +05:30
! == cleaned:end ===================================================================================
2020-12-20 21:41:43 +05:30
module function phase_thermal_constitutive(Delta_t,ph,en) result(status)
real(pREAL), intent(in) :: Delta_t
2021-04-25 11:36:52 +05:30
integer, intent(in) :: ph, en
integer(kind(STATUS_OK)) :: status
2021-07-17 00:20:08 +05:30
end function phase_thermal_constitutive
2021-01-08 02:45:18 +05:30
module function phase_damage_constitutive(Delta_t,co,ce) result(status)
real(pREAL), intent(in) :: Delta_t
2022-02-04 16:55:11 +05:30
integer, intent(in) :: co, ce
integer(kind(STATUS_OK)) :: status
end function phase_damage_constitutive
2021-01-08 02:45:18 +05:30
module function phase_mechanical_constitutive(Delta_t,co,ce) result(status)
real(pREAL), intent(in) :: Delta_t
2022-02-04 16:55:11 +05:30
integer, intent(in) :: co, ce
integer(kind(STATUS_OK)) :: status
2021-07-17 00:20:08 +05:30
end function phase_mechanical_constitutive
2020-12-24 14:52:41 +05:30
2021-07-17 15:20:21 +05:30
!ToDo: Merge all the stiffness functions
2021-11-18 21:07:34 +05:30
module function phase_homogenizedC66(ph,en) result(C)
2021-04-13 00:51:15 +05:30
integer, intent(in) :: ph, en
real(pREAL), dimension(6,6) :: C
2021-11-18 21:07:34 +05:30
end function phase_homogenizedC66
2021-11-19 02:29:09 +05:30
module function phase_damage_C66(C66,ph,en) result(C66_degraded)
real(pREAL), dimension(6,6), intent(in) :: C66
2021-11-18 21:26:36 +05:30
integer, intent(in) :: ph,en
real(pREAL), dimension(6,6) :: C66_degraded
2021-11-18 21:26:36 +05:30
end function phase_damage_C66
2021-04-09 03:10:20 +05:30
module function phase_f_phi(phi,co,ce) result(f)
2021-04-07 16:32:42 +05:30
integer, intent(in) :: ce,co
real(pREAL), intent(in) :: &
2020-07-15 18:05:21 +05:30
phi !< damage parameter
real(pREAL) :: &
2021-04-09 03:10:20 +05:30
f
2021-04-08 17:01:21 +05:30
end function phase_f_phi
2020-07-15 18:05:21 +05:30
2021-04-25 11:36:52 +05:30
module function phase_f_T(ph,en) result(f)
integer, intent(in) :: ph, en
real(pREAL) :: f
2021-04-08 02:11:49 +05:30
end function phase_f_T
2020-07-15 18:05:21 +05:30
2022-02-04 22:12:05 +05:30
module subroutine plastic_dependentState(ph,en)
2020-07-15 18:05:21 +05:30
integer, intent(in) :: &
2022-02-04 22:12:05 +05:30
ph, &
en
2021-01-26 16:47:00 +05:30
end subroutine plastic_dependentState
2020-07-15 18:05:21 +05:30
2021-02-13 23:27:41 +05:30
module subroutine damage_anisobrittle_LiAndItsTangent(L_i, dL_i_dM_i, M_i, ph,en)
integer, intent(in) :: ph, en
real(pREAL), intent(in), dimension(3,3) :: &
M_i
real(pREAL), intent(out), dimension(3,3) :: &
L_i !< damage velocity gradient
real(pREAL), intent(out), dimension(3,3,3,3) :: &
dL_i_dM_i !< derivative of L_i with respect to M_i
2021-04-06 15:48:48 +05:30
end subroutine damage_anisobrittle_LiAndItsTangent
2021-02-13 23:27:41 +05:30
2020-12-21 12:35:38 +05:30
end interface
2019-06-15 19:10:22 +05:30
public :: &
2021-02-09 03:51:53 +05:30
phase_init, &
2021-11-18 21:07:34 +05:30
phase_homogenizedC66, &
2021-04-08 17:01:21 +05:30
phase_f_phi, &
2021-04-08 02:11:49 +05:30
phase_f_T, &
2021-04-11 16:51:27 +05:30
phase_K_phi, &
phase_K_T, &
phase_mu_phi, &
phase_mu_T, &
2023-01-19 22:07:45 +05:30
phase_result, &
2021-02-09 03:51:53 +05:30
phase_allocateState, &
phase_forward, &
phase_restore, &
2020-12-21 15:27:18 +05:30
converged, &
2021-07-17 00:20:08 +05:30
phase_mechanical_constitutive, &
phase_thermal_constitutive, &
phase_damage_constitutive, &
2021-02-09 03:51:53 +05:30
phase_mechanical_dPdF, &
crystallite_orientations, &
crystallite_push33ToRef, &
2021-02-09 03:51:53 +05:30
phase_restartWrite, &
phase_restartRead, &
phase_thermal_setField, &
2021-04-08 17:01:21 +05:30
phase_set_phi, &
2021-04-08 00:36:29 +05:30
phase_P, &
phase_set_F, &
phase_F
2020-12-21 15:27:18 +05:30
contains
!--------------------------------------------------------------------------------------------------
2021-04-06 13:44:52 +05:30
!> @brief Initialize constitutive models for individual physics
!--------------------------------------------------------------------------------------------------
2023-12-28 23:48:25 +05:30
subroutine phase_init()
2019-06-15 19:10:22 +05:30
integer :: &
2021-05-22 20:51:07 +05:30
ph, ce, co, ma
type(tDict), pointer :: &
phases, &
phase, &
num_phase, &
num_mech
character(len=:), allocatable :: refs
2020-07-02 02:21:21 +05:30
2020-12-23 19:33:03 +05:30
print'(/,1x,a)', '<<<+- phase init -+>>>'; flush(IO_STDOUT)
phases => config_material%get_dict('phase')
allocate(phase_lattice(phases%length))
allocate(phase_cOverA(phases%length),source=-1.0_pREAL)
allocate(phase_rho(phases%length))
2021-07-24 18:33:26 +05:30
allocate(phase_O_0(phases%length))
2021-05-22 20:51:07 +05:30
do ph = 1,phases%length
print'(/,1x,a,i0,a)', 'phase ',ph,': '//phases%key(ph)
phase => phases%get_dict(ph)
refs = config_listReferences(phase,indent=3)
if (len(refs) > 0) print'(/,1x,a)', refs
phase_rho(ph) = phase%get_asReal('rho',defaultVal=0.0_pREAL)
2023-06-04 10:47:38 +05:30
phase_lattice(ph) = phase%get_asStr('lattice')
if (all(phase_lattice(ph) /= ['cF','cI','hP','tI'])) &
2023-06-04 10:47:38 +05:30
call IO_error(130,ext_msg='phase_init: '//phase%get_asStr('lattice'))
if (any(phase_lattice(ph) == ['hP','tI'])) &
phase_cOverA(ph) = phase%get_asReal('c/a')
2023-01-23 13:01:59 +05:30
allocate(phase_O_0(ph)%data(count(material_ID_phase==ph)))
end do
2021-05-22 20:51:07 +05:30
2023-01-23 13:01:59 +05:30
do ce = 1, size(material_ID_phase,2)
2021-05-22 20:51:07 +05:30
ma = discretization_materialAt((ce-1)/discretization_nIPs+1)
2023-01-23 13:01:59 +05:30
do co = 1,homogenization_Nconstituents(material_ID_homogenization(ce))
ph = material_ID_phase(co,ce)
phase_O_0(ph)%data(material_entry_phase(co,ce)) = material_O_0(ma)%data(co)
end do
end do
2021-05-22 20:51:07 +05:30
2021-07-24 18:33:26 +05:30
allocate(phase_O(phases%length))
2021-05-22 20:51:07 +05:30
do ph = 1,phases%length
2021-07-24 18:33:26 +05:30
phase_O(ph)%data = phase_O_0(ph)%data
end do
2021-05-22 20:51:07 +05:30
num_phase => config_numerics%get_dict('phase',defaultVal=emptyDict)
num_mech => num_phase%get_dict('mechanical', defaultVal=emptyDict)
call mechanical_init(phases,num_mech)
call damage_init()
call thermal_init(phases)
call crystallite_init()
2021-02-09 03:51:53 +05:30
end subroutine phase_init
2009-03-06 15:32:36 +05:30
2020-08-15 19:32:10 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Allocate the components of the state structure for a given phase
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_allocateState(state, &
2022-02-04 03:27:32 +05:30
NEntries,sizeState,sizeDotState,sizeDeltaState,offsetDeltaState)
2020-08-15 19:32:10 +05:30
class(tState), intent(inout) :: &
2020-08-15 19:32:10 +05:30
state
integer, intent(in) :: &
2021-04-06 13:44:52 +05:30
NEntries, &
2020-08-15 19:32:10 +05:30
sizeState, &
sizeDotState, &
sizeDeltaState
2022-02-04 03:27:32 +05:30
integer, intent(in), optional :: &
offsetDeltaState
2020-12-29 02:11:48 +05:30
2020-08-15 19:32:10 +05:30
state%sizeState = sizeState
state%sizeDotState = sizeDotState
state%sizeDeltaState = sizeDeltaState
2022-02-04 03:27:32 +05:30
if (present(offsetDeltaState)) then
state%offsetDeltaState = offsetDeltaState ! ToDo: this is a fix for broken nonlocal
else
state%offsetDeltaState = sizeState-sizeDeltaState ! deltaState occupies latter part of state by definition
end if
2020-08-15 19:32:10 +05:30
allocate(state%atol (sizeState), source=0.0_pREAL)
allocate(state%state0 (sizeState,NEntries), source=0.0_pREAL)
allocate(state%state (sizeState,NEntries), source=0.0_pREAL)
2020-08-15 19:32:10 +05:30
allocate(state%dotState (sizeDotState,NEntries), source=0.0_pREAL)
2020-08-15 19:32:10 +05:30
allocate(state%deltaState (sizeDeltaState,NEntries), source=0.0_pREAL)
state%deltaState2 => state%state(state%offsetDeltaState+1: &
state%offsetDeltaState+state%sizeDeltaState,:)
2020-08-15 19:32:10 +05:30
2021-02-09 03:51:53 +05:30
end subroutine phase_allocateState
2020-08-15 19:32:10 +05:30
2020-12-20 13:57:37 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Restore data after homog cutback.
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_restore(ce,includeL)
2020-12-20 13:57:37 +05:30
2020-12-28 02:15:31 +05:30
logical, intent(in) :: includeL
2021-01-24 14:09:40 +05:30
integer, intent(in) :: ce
2020-12-28 02:15:31 +05:30
2020-12-20 13:57:37 +05:30
2021-02-09 03:51:53 +05:30
call mechanical_restore(ce,includeL)
2021-07-17 03:02:08 +05:30
call damage_restore(ce)
2020-12-28 02:15:31 +05:30
2021-02-09 03:51:53 +05:30
end subroutine phase_restore
2020-12-20 13:57:37 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_forward()
2020-12-20 13:57:37 +05:30
2021-02-09 03:51:53 +05:30
call mechanical_forward()
2021-04-06 13:52:47 +05:30
call damage_forward()
2020-12-31 14:24:13 +05:30
call thermal_forward()
2020-12-23 11:28:54 +05:30
2021-02-09 03:51:53 +05:30
end subroutine phase_forward
2020-12-20 13:57:37 +05:30
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes constitutive results to HDF5 output file
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
2023-01-19 22:07:45 +05:30
subroutine phase_result()
2019-12-19 00:35:51 +05:30
2020-12-22 13:15:01 +05:30
integer :: ph
character(len=:), allocatable :: group
2023-01-19 22:07:45 +05:30
call result_closeGroup(result_addGroup('/current/phase/'))
2020-12-22 13:15:01 +05:30
do ph = 1, size(material_name_phase)
2020-12-22 16:50:00 +05:30
group = '/current/phase/'//trim(material_name_phase(ph))//'/'
2023-01-19 22:07:45 +05:30
call result_closeGroup(result_addGroup(group))
2020-12-22 13:15:01 +05:30
2023-01-19 22:07:45 +05:30
call mechanical_result(group,ph)
call damage_result(group,ph)
call thermal_result(group,ph)
2020-12-22 13:15:01 +05:30
end do
2020-12-22 13:15:01 +05:30
2023-01-19 22:07:45 +05:30
end subroutine phase_result
2018-12-05 04:25:39 +05:30
2020-08-15 19:32:10 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Allocate and initialize.
!--------------------------------------------------------------------------------------------------
2020-12-31 14:24:13 +05:30
subroutine crystallite_init()
integer :: &
2021-05-23 13:40:25 +05:30
ce, &
2020-12-29 02:11:48 +05:30
co, & !< counter in integration point component loop
2022-02-04 13:19:00 +05:30
en, ph
2023-12-29 00:21:22 +05:30
!$OMP PARALLEL DO PRIVATE(ph,en)
do ce = 1, size(material_ID_homogenization)
do co = 1,homogenization_Nconstituents(material_ID_homogenization(ce))
ph = material_ID_phase(co,ce)
en = material_entry_phase(co,ce)
call crystallite_orientations(co,ce)
call plastic_dependentState(ph,en) ! update dependent state variables to be consistent with basic states
end do
end do
!$OMP END PARALLEL DO
end subroutine crystallite_init
!--------------------------------------------------------------------------------------------------
2023-12-29 00:21:22 +05:30
!> @brief Update orientations and, if needed, compatibility.
!--------------------------------------------------------------------------------------------------
2023-12-29 00:21:22 +05:30
subroutine crystallite_orientations(co,ce)
2020-12-27 16:18:02 +05:30
integer, intent(in) :: &
2023-12-29 00:21:22 +05:30
co, &
ce
2020-12-23 11:44:07 +05:30
2021-05-22 22:12:18 +05:30
integer :: ph, en
2023-12-29 00:21:22 +05:30
ph = material_ID_phase(co,ce)
en = material_entry_phase(co,ce)
2021-07-24 18:33:26 +05:30
call phase_O(ph)%data(en)%fromMatrix(transpose(math_rotationalPart(mechanical_F_e(ph,en))))
2020-12-27 16:18:02 +05:30
end subroutine crystallite_orientations
!--------------------------------------------------------------------------------------------------
2023-12-29 00:21:22 +05:30
!> @brief Map 2nd order tensor to reference configuration.
!--------------------------------------------------------------------------------------------------
function crystallite_push33ToRef(co,ce, tensor33)
real(pREAL), dimension(3,3), intent(in) :: tensor33
integer, intent(in):: &
co, &
ce
real(pREAL), dimension(3,3) :: crystallite_push33ToRef
2020-12-30 04:44:48 +05:30
real(pREAL), dimension(3,3) :: T
2021-04-06 15:08:44 +05:30
integer :: ph, en
2022-06-24 10:34:52 +05:30
2023-01-23 13:01:59 +05:30
ph = material_ID_phase(co,ce)
en = material_entry_phase(co,ce)
2021-07-24 18:33:26 +05:30
T = matmul(phase_O_0(ph)%data(en)%asMatrix(),transpose(math_inv33(phase_F(co,ce)))) ! ToDo: initial orientation correct?
2020-12-30 14:24:06 +05:30
crystallite_push33ToRef = matmul(transpose(T),matmul(tensor33,T))
end function crystallite_push33ToRef
!--------------------------------------------------------------------------------------------------
2023-12-29 01:58:57 +05:30
!> @brief Determine whether a point is converged.
!--------------------------------------------------------------------------------------------------
logical pure function converged(residuum,state,atol)
2023-12-29 00:21:22 +05:30
real(pREAL), intent(in), dimension(:) :: &
residuum, state, atol
2023-12-29 00:21:22 +05:30
real(pREAL) :: &
rTol
2023-12-29 00:21:22 +05:30
rTol = num%rTol_crystalliteState
converged = all(abs(residuum) <= max(atol, rtol*abs(state)))
end function converged
!--------------------------------------------------------------------------------------------------
2021-04-06 13:44:52 +05:30
!> @brief Write restart data to file.
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_restartWrite(fileHandle)
2020-12-30 02:01:22 +05:30
integer(HID_T), intent(in) :: fileHandle
2020-12-30 04:44:48 +05:30
integer(HID_T), dimension(2) :: groupHandle
2020-12-28 13:37:35 +05:30
integer :: ph
2020-12-30 04:44:48 +05:30
groupHandle(1) = HDF5_addGroup(fileHandle,'phase')
do ph = 1, size(material_name_phase)
groupHandle(2) = HDF5_addGroup(groupHandle(1),material_name_phase(ph))
2021-02-09 03:51:53 +05:30
call mechanical_restartWrite(groupHandle(2),ph)
2022-01-11 01:04:50 +05:30
call thermal_restartWrite(groupHandle(2),ph)
call damage_restartWrite(groupHandle(2),ph)
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
end do
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(1))
2021-02-09 03:51:53 +05:30
end subroutine phase_restartWrite
!--------------------------------------------------------------------------------------------------
2021-04-06 13:44:52 +05:30
!> @brief Read restart data from file.
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_restartRead(fileHandle)
2020-12-30 02:01:22 +05:30
integer(HID_T), intent(in) :: fileHandle
2020-12-30 04:44:48 +05:30
integer(HID_T), dimension(2) :: groupHandle
2020-12-28 13:37:35 +05:30
integer :: ph
2020-12-30 04:44:48 +05:30
groupHandle(1) = HDF5_openGroup(fileHandle,'phase')
do ph = 1, size(material_name_phase)
groupHandle(2) = HDF5_openGroup(groupHandle(1),material_name_phase(ph))
2021-02-09 03:51:53 +05:30
call mechanical_restartRead(groupHandle(2),ph)
2022-01-11 01:04:50 +05:30
call thermal_restartRead(groupHandle(2),ph)
call damage_restartRead(groupHandle(2),ph)
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
end do
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(1))
2021-02-09 03:51:53 +05:30
end subroutine phase_restartRead
end module phase