DAMASK_EICMD/src/phase.f90

686 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
2019-06-15 19:10:22 +05:30
use math
use rotations
2019-06-15 19:10:22 +05:30
use IO
use config
use material
use results
use lattice
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
2019-06-15 19:10:22 +05:30
implicit none
private
2020-12-23 22:00:19 +05:30
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 :: &
phase_orientation0, &
phase_orientation
2020-12-21 00:50:39 +05:30
type :: tTensorContainer
real(pReal), dimension(:,:,:), allocatable :: data
end type
type :: tNumerics
integer :: &
iJacoLpresiduum, & !< frequency of Jacobian update of residuum in Lp
nState, & !< state loop limit
nStress !< stress loop limit
real(pReal) :: &
subStepMinCryst, & !< minimum (relative) size of sub-step allowed during cutback
subStepSizeCryst, & !< size of first substep when cutback
subStepSizeLp, & !< size of first substep when cutback in Lp calculation
subStepSizeLi, & !< size of first substep when cutback in Li calculation
stepIncreaseCryst, & !< increase of next substep size when previous substep converged
rtol_crystalliteState, & !< relative tolerance in state loop
rtol_crystalliteStress, & !< relative tolerance in stress loop
atol_crystalliteStress !< absolute tolerance in stress loop
end type tNumerics
type(tNumerics) :: num ! numerics parameters. Better name?
type :: tDebugOptions
logical :: &
basic, &
extensive, &
selective
integer :: &
element, &
ip, &
grain
end type tDebugOptions
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(materials,phases)
class(tNode), pointer :: materials,phases
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)
class(tNode), pointer :: phases
end subroutine thermal_init
2020-12-21 16:44:09 +05:30
2021-02-09 03:51:53 +05:30
module subroutine mechanical_results(group,ph)
2020-12-21 16:44:09 +05:30
character(len=*), intent(in) :: group
integer, intent(in) :: ph
2021-02-09 03:51:53 +05:30
end subroutine mechanical_results
2020-12-21 16:44:09 +05:30
2020-12-22 16:50:00 +05:30
module subroutine damage_results(group,ph)
character(len=*), intent(in) :: group
integer, intent(in) :: ph
end subroutine damage_results
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
2020-12-30 02:01:22 +05:30
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
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
2020-12-30 14:24:06 +05:30
2021-04-13 00:51:15 +05:30
module function mechanical_S(ph,en) result(S)
integer, intent(in) :: ph,en
2020-12-30 14:24:06 +05:30
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
2020-12-30 17:15:08 +05:30
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
2020-12-30 14:24:06 +05:30
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
2021-04-08 17:01:21 +05:30
module function phase_F(co,ce) result(F)
integer, intent(in) :: co, ce
2021-04-08 17:01:21 +05:30
real(pReal), dimension(3,3) :: F
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
2020-12-30 15:33:13 +05:30
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
2021-04-25 11:36:52 +05:30
module function thermal_T(ph,en) result(T)
integer, intent(in) :: ph,en
2020-12-30 14:24:06 +05:30
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
2021-02-13 11:45:57 +05:30
real(pReal) :: phi
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)
2020-12-30 14:24:06 +05:30
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)
2021-01-24 17:56:01 +05:30
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)
2021-01-21 01:24:31 +05:30
real(pReal), intent(in) :: phi
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
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
end function phase_K_phi
module function phase_mu_T(co,ce) result(mu)
integer, intent(in) :: co, ce
real(pReal) :: mu
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
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
2021-07-17 00:20:08 +05:30
module function phase_thermal_constitutive(Delta_t,ph,en) result(converged_)
real(pReal), intent(in) :: Delta_t
2021-04-25 11:36:52 +05:30
integer, intent(in) :: ph, en
logical :: converged_
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,ip,el) result(converged_)
2021-07-17 15:20:21 +05:30
real(pReal), intent(in) :: Delta_t
integer, intent(in) :: co, ip, el
logical :: converged_
end function phase_damage_constitutive
2021-01-08 02:45:18 +05:30
2021-07-17 15:20:21 +05:30
module function phase_mechanical_constitutive(Delta_t,co,ip,el) result(converged_)
real(pReal), intent(in) :: Delta_t
2020-12-24 14:52:41 +05:30
integer, intent(in) :: co, ip, el
2020-12-27 20:43:31 +05:30
logical :: converged_
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-04-13 00:51:15 +05:30
module function phase_homogenizedC(ph,en) result(C)
integer, intent(in) :: ph, en
2020-12-23 22:00:19 +05:30
real(pReal), dimension(6,6) :: C
2021-02-09 03:51:53 +05:30
end function phase_homogenizedC
2021-05-22 01:10:46 +05:30
module function phase_damage_C(C_homogenized,ph,en) result(C)
real(pReal), dimension(3,3,3,3), intent(in) :: C_homogenized
integer, intent(in) :: ph,en
real(pReal), dimension(3,3,3,3) :: C
end function phase_damage_C
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
2020-07-15 18:05:21 +05:30
real(pReal), intent(in) :: &
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
2021-04-09 03:10:20 +05:30
real(pReal) :: f
2021-04-08 02:11:49 +05:30
end function phase_f_T
2020-07-15 18:05:21 +05:30
module subroutine plastic_nonlocal_updateCompatibility(orientation,ph,i,e)
2020-07-10 21:59:36 +05:30
integer, intent(in) :: &
ph, &
2020-07-10 21:59:36 +05:30
i, &
e
2021-05-22 22:12:18 +05:30
type(tRotationContainer), dimension(:), intent(in) :: orientation
2020-07-10 21:59:36 +05:30
end subroutine plastic_nonlocal_updateCompatibility
2021-01-26 16:47:00 +05:30
module subroutine plastic_dependentState(co,ip,el)
2020-07-15 18:05:21 +05:30
integer, intent(in) :: &
2021-01-08 12:07:51 +05:30
co, & !< component-ID of integration point
2020-07-15 18:05:21 +05:30
ip, & !< integration point
el !< element
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(Ld, dLd_dTstar, S, ph,en)
integer, intent(in) :: ph, en
2021-02-13 23:27:41 +05:30
real(pReal), intent(in), dimension(3,3) :: &
S
real(pReal), intent(out), dimension(3,3) :: &
Ld !< damage velocity gradient
real(pReal), intent(out), dimension(3,3,3,3) :: &
dLd_dTstar !< derivative of Ld with respect to Tstar (4th-order tensor)
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
2020-07-02 02:21:21 +05:30
type(tDebugOptions) :: debugConstitutive
#if __INTEL_COMPILER >= 1900
public :: &
prec, &
math, &
rotations, &
IO, &
config, &
material, &
results, &
lattice, &
discretization, &
HDF5_utilities
#endif
2019-06-15 19:10:22 +05:30
public :: &
2021-02-09 03:51:53 +05:30
phase_init, &
phase_homogenizedC, &
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, &
2021-02-09 03:51:53 +05:30
phase_results, &
phase_allocateState, &
phase_forward, &
phase_restore, &
2020-08-15 19:32:10 +05:30
plastic_nonlocal_updateCompatibility, &
2020-12-21 15:27:18 +05:30
converged, &
crystallite_init, &
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
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +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
2020-07-02 02:21:21 +05:30
class (tNode), pointer :: &
2020-08-15 19:32:10 +05:30
debug_constitutive, &
materials, &
phases, &
phase
2020-07-02 02:21:21 +05:30
2020-12-23 19:33:03 +05:30
2021-01-27 15:14:03 +05:30
print'(/,a)', ' <<<+- phase init -+>>>'; flush(IO_STDOUT)
2021-03-26 22:07:38 +05:30
debug_constitutive => config_debug%get('phase', defaultVal=emptyList)
debugConstitutive%basic = debug_constitutive%contains('basic')
debugConstitutive%extensive = debug_constitutive%contains('extensive')
debugConstitutive%selective = debug_constitutive%contains('selective')
debugConstitutive%element = config_debug%get_asInt('element', defaultVal = 1)
debugConstitutive%ip = config_debug%get_asInt('integrationpoint',defaultVal = 1)
debugConstitutive%grain = config_debug%get_asInt('constituent', defaultVal = 1)
2020-07-02 02:21:21 +05:30
materials => config_material%get('material')
2021-03-26 22:07:38 +05:30
phases => config_material%get('phase')
allocate(phase_lattice(phases%length))
allocate(phase_cOverA(phases%length),source=-1.0_pReal)
allocate(phase_rho(phases%length))
2021-05-22 20:51:07 +05:30
allocate(phase_orientation0(phases%length))
2021-05-22 20:51:07 +05:30
do ph = 1,phases%length
phase => phases%get(ph)
phase_lattice(ph) = phase%get_asString('lattice')
if (all(phase_lattice(ph) /= ['cF','cI','hP','tI'])) &
call IO_error(130,ext_msg='phase_init: '//phase%get_asString('lattice'))
if (any(phase_lattice(ph) == ['hP','tI'])) &
phase_cOverA(ph) = phase%get_asFloat('c/a')
phase_rho(ph) = phase%get_asFloat('rho',defaultVal=0.0_pReal)
2021-05-22 20:51:07 +05:30
allocate(phase_orientation0(ph)%data(count(material_phaseID==ph)))
enddo
do ce = 1, size(material_phaseID,2)
ma = discretization_materialAt((ce-1)/discretization_nIPs+1)
do co = 1,homogenization_Nconstituents(material_homogenizationID(ce))
ph = material_phaseID(co,ce)
2021-05-22 22:12:18 +05:30
phase_orientation0(ph)%data(material_phaseEntry(co,ce)) = material_orientation0(ma)%data(co)
2021-05-22 20:51:07 +05:30
enddo
enddo
allocate(phase_orientation(phases%length))
do ph = 1,phases%length
phase_orientation(ph)%data = phase_orientation0(ph)%data
enddo
call mechanical_init(materials,phases)
call damage_init
call thermal_init(phases)
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, &
2021-04-06 13:44:52 +05:30
NEntries,sizeState,sizeDotState,sizeDeltaState)
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
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
state%offsetDeltaState = sizeState-sizeDeltaState ! deltaState occupies latter part of state by definition
2021-04-06 13:44:52 +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
2021-04-06 13:44:52 +05:30
allocate(state%dotState (sizeDotState,NEntries), source=0.0_pReal)
2020-08-15 19:32:10 +05:30
2021-04-06 13:44:52 +05:30
allocate(state%deltaState (sizeDeltaState,NEntries), source=0.0_pReal)
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
!--------------------------------------------------------------------------------------------------
2021-02-09 03:51:53 +05:30
subroutine phase_results()
2019-12-19 00:35:51 +05:30
2020-12-22 13:15:01 +05:30
integer :: ph
character(len=:), allocatable :: group
2020-12-22 16:50:00 +05:30
call results_closeGroup(results_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))//'/'
2020-12-22 13:15:01 +05:30
call results_closeGroup(results_addGroup(group))
2021-02-09 03:51:53 +05:30
call mechanical_results(group,ph)
2020-12-22 16:50:00 +05:30
call damage_results(group,ph)
2020-12-22 13:15:01 +05:30
enddo
2021-02-09 03:51:53 +05:30
end subroutine phase_results
2018-12-05 04:25:39 +05:30
2020-08-15 19:32:10 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief allocates and initialize per grain variables
!--------------------------------------------------------------------------------------------------
2020-12-31 14:24:13 +05:30
subroutine crystallite_init()
integer :: &
2020-12-23 17:52:11 +05:30
ph, &
2021-05-23 13:40:25 +05:30
ce, &
2020-12-29 02:11:48 +05:30
co, & !< counter in integration point component loop
ip, & !< counter in integration point loop
el, & !< counter in element loop
cMax, & !< maximum number of integration point components
iMax, & !< maximum number of integration points
eMax !< maximum number of elements
2021-01-08 02:45:18 +05:30
class(tNode), pointer :: &
num_crystallite, &
phases
2020-12-23 11:44:07 +05:30
print'(/,a)', ' <<<+- crystallite init -+>>>'
cMax = homogenization_maxNconstituents
iMax = discretization_nIPs
eMax = discretization_Nelems
num_crystallite => config_numerics%get('crystallite',defaultVal=emptyDict)
num%subStepMinCryst = num_crystallite%get_asFloat ('subStepMin', defaultVal=1.0e-3_pReal)
num%subStepSizeCryst = num_crystallite%get_asFloat ('subStepSize', defaultVal=0.25_pReal)
num%stepIncreaseCryst = num_crystallite%get_asFloat ('stepIncrease', defaultVal=1.5_pReal)
num%subStepSizeLp = num_crystallite%get_asFloat ('subStepSizeLp', defaultVal=0.5_pReal)
num%subStepSizeLi = num_crystallite%get_asFloat ('subStepSizeLi', defaultVal=0.5_pReal)
num%rtol_crystalliteState = num_crystallite%get_asFloat ('rtol_State', defaultVal=1.0e-6_pReal)
num%rtol_crystalliteStress = num_crystallite%get_asFloat ('rtol_Stress', defaultVal=1.0e-6_pReal)
num%atol_crystalliteStress = num_crystallite%get_asFloat ('atol_Stress', defaultVal=1.0e-8_pReal)
num%iJacoLpresiduum = num_crystallite%get_asInt ('iJacoLpresiduum', defaultVal=1)
num%nState = num_crystallite%get_asInt ('nState', defaultVal=20)
num%nStress = num_crystallite%get_asInt ('nStress', defaultVal=40)
if(num%subStepMinCryst <= 0.0_pReal) call IO_error(301,ext_msg='subStepMinCryst')
if(num%subStepSizeCryst <= 0.0_pReal) call IO_error(301,ext_msg='subStepSizeCryst')
if(num%stepIncreaseCryst <= 0.0_pReal) call IO_error(301,ext_msg='stepIncreaseCryst')
if(num%subStepSizeLp <= 0.0_pReal) call IO_error(301,ext_msg='subStepSizeLp')
if(num%subStepSizeLi <= 0.0_pReal) call IO_error(301,ext_msg='subStepSizeLi')
if(num%rtol_crystalliteState <= 0.0_pReal) call IO_error(301,ext_msg='rtol_crystalliteState')
if(num%rtol_crystalliteStress <= 0.0_pReal) call IO_error(301,ext_msg='rtol_crystalliteStress')
if(num%atol_crystalliteStress <= 0.0_pReal) call IO_error(301,ext_msg='atol_crystalliteStress')
if(num%iJacoLpresiduum < 1) call IO_error(301,ext_msg='iJacoLpresiduum')
if(num%nState < 1) call IO_error(301,ext_msg='nState')
if(num%nStress< 1) call IO_error(301,ext_msg='nStress')
phases => config_material%get('phase')
2020-12-20 22:52:04 +05:30
print'(a42,1x,i10)', ' # of elements: ', eMax
print'(a42,1x,i10)', ' # of integration points/element: ', iMax
print'(a42,1x,i10)', 'max # of constituents/integration point: ', cMax
flush(IO_STDOUT)
2021-05-23 22:14:24 +05:30
!$OMP PARALLEL DO PRIVATE(ce)
2021-05-23 03:40:46 +05:30
do el = 1, eMax
do ip = 1, iMax
2021-05-23 13:40:25 +05:30
ce = (el-1)*discretization_nIPs + ip
do co = 1,homogenization_Nconstituents(material_homogenizationID(ce))
2020-12-27 16:18:02 +05:30
call crystallite_orientations(co,ip,el)
call plastic_dependentState(co,ip,el) ! update dependent state variables to be consistent with basic states
enddo
enddo
enddo
!$OMP END PARALLEL DO
end subroutine crystallite_init
!--------------------------------------------------------------------------------------------------
!> @brief calculates orientations
!--------------------------------------------------------------------------------------------------
2020-12-27 16:18:02 +05:30
subroutine crystallite_orientations(co,ip,el)
2020-12-27 16:18:02 +05:30
integer, intent(in) :: &
co, & !< counter in integration point component loop
ip, & !< counter in integration point loop
el !< counter in element loop
2020-12-23 11:44:07 +05:30
2021-05-22 22:12:18 +05:30
integer :: ph, en
ph = material_phaseID(co,(el-1)*discretization_nIPs + ip)
en = material_phaseEntry(co,(el-1)*discretization_nIPs + ip)
2021-05-22 22:12:18 +05:30
call phase_orientation(ph)%data(en)%fromMatrix(transpose(math_rotationalPart(mechanical_F_e(ph,en))))
2020-12-27 16:18:02 +05:30
if (plasticState(material_phaseAt(1,el))%nonlocal) &
2021-05-22 22:12:18 +05:30
call plastic_nonlocal_updateCompatibility(phase_orientation,material_phaseAt(1,el),ip,el)
end subroutine crystallite_orientations
!--------------------------------------------------------------------------------------------------
!> @brief Map 2nd order tensor to reference config
!--------------------------------------------------------------------------------------------------
function crystallite_push33ToRef(co,ce, tensor33)
real(pReal), dimension(3,3), intent(in) :: tensor33
integer, intent(in):: &
co, &
ce
2020-12-23 12:42:56 +05:30
real(pReal), dimension(3,3) :: crystallite_push33ToRef
2020-12-30 04:44:48 +05:30
2020-12-29 22:57:24 +05:30
real(pReal), dimension(3,3) :: T
2021-04-06 15:08:44 +05:30
integer :: ph, en
2021-04-06 15:08:44 +05:30
ph = material_phaseID(co,ce)
en = material_phaseEntry(co,ce)
2021-05-22 20:51:07 +05:30
T = matmul(phase_orientation0(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
!--------------------------------------------------------------------------------------------------
!> @brief determines whether a point is converged
!--------------------------------------------------------------------------------------------------
logical pure function converged(residuum,state,atol)
real(pReal), intent(in), dimension(:) ::&
residuum, state, atol
real(pReal) :: &
rTol
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)
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
enddo
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)
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
enddo
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