DAMASK_EICMD/src/constitutive_thermal.f90

169 lines
5.8 KiB
Fortran
Raw Normal View History

2020-07-15 18:05:21 +05:30
!----------------------------------------------------------------------------------------------------
2020-08-13 00:44:00 +05:30
!> @brief internal microstructure state for all thermal sources and kinematics constitutive models
2020-07-15 18:05:21 +05:30
!----------------------------------------------------------------------------------------------------
submodule(constitutive) constitutive_thermal
type :: tDataContainer
real(pReal), dimension(:), allocatable :: T
end type tDataContainer
type(tDataContainer), dimension(:), allocatable :: current
interface
2020-08-15 19:32:10 +05:30
module function source_thermal_dissipation_init(source_length) result(mySources)
integer, intent(in) :: source_length
logical, dimension(:,:), allocatable :: mySources
end function source_thermal_dissipation_init
2020-12-30 14:24:06 +05:30
2020-08-15 19:32:10 +05:30
module function source_thermal_externalheat_init(source_length) result(mySources)
integer, intent(in) :: source_length
logical, dimension(:,:), allocatable :: mySources
end function source_thermal_externalheat_init
2020-08-13 00:44:00 +05:30
2020-08-15 19:32:10 +05:30
module function kinematics_thermal_expansion_init(kinematics_length) result(myKinematics)
integer, intent(in) :: kinematics_length
logical, dimension(:,:), allocatable :: myKinematics
end function kinematics_thermal_expansion_init
2020-07-10 18:29:07 +05:30
module subroutine source_thermal_dissipation_getRateAndItsTangent(TDot, dTDot_dT, Tstar, Lp, phase)
integer, intent(in) :: &
2020-07-15 18:05:21 +05:30
phase !< phase ID of element
real(pReal), intent(in), dimension(3,3) :: &
2020-08-13 00:44:00 +05:30
Tstar !< 2nd Piola Kirchhoff stress tensor for a given element
real(pReal), intent(in), dimension(3,3) :: &
2020-07-15 18:05:21 +05:30
Lp !< plastic velocuty gradient for a given element
real(pReal), intent(out) :: &
TDot, &
dTDot_dT
end subroutine source_thermal_dissipation_getRateAndItsTangent
module subroutine source_thermal_externalheat_getRateAndItsTangent(TDot, dTDot_dT, phase, of)
integer, intent(in) :: &
phase, &
of
real(pReal), intent(out) :: &
TDot, &
dTDot_dT
end subroutine source_thermal_externalheat_getRateAndItsTangent
end interface
2020-07-12 20:14:26 +05:30
contains
2020-07-12 20:14:26 +05:30
!----------------------------------------------------------------------------------------------
!< @brief initializes thermal sources and kinematics mechanism
!----------------------------------------------------------------------------------------------
module subroutine thermal_init(phases)
class(tNode), pointer :: &
phases
integer :: &
ph, &
Nconstituents
print'(/,a)', ' <<<+- constitutive_mech init -+>>>'
allocate(current(phases%length))
do ph = 1, phases%length
Nconstituents = count(material_phaseAt == ph) * discretization_nIPs
allocate(current(ph)%T(Nconstituents))
enddo
! initialize source mechanisms
2020-08-15 19:32:10 +05:30
if(maxval(phase_Nsources) /= 0) then
where(source_thermal_dissipation_init (maxval(phase_Nsources))) phase_source = SOURCE_thermal_dissipation_ID
where(source_thermal_externalheat_init(maxval(phase_Nsources))) phase_source = SOURCE_thermal_externalheat_ID
2020-12-30 14:24:06 +05:30
endif
!--------------------------------------------------------------------------------------------------
!initialize kinematic mechanisms
2020-08-15 19:32:10 +05:30
if(maxval(phase_Nkinematics) /= 0) where(kinematics_thermal_expansion_init(maxval(phase_Nkinematics))) &
phase_kinematics = KINEMATICS_thermal_expansion_ID
end subroutine thermal_init
2020-07-12 20:14:26 +05:30
!----------------------------------------------------------------------------------------------
!< @brief calculates thermal dissipation rate
!----------------------------------------------------------------------------------------------
2020-12-29 11:50:37 +05:30
module subroutine constitutive_thermal_getRateAndItsTangents(TDot, dTDot_dT, T, ip, el)
2020-07-12 20:14:26 +05:30
integer, intent(in) :: &
2020-07-15 18:05:21 +05:30
ip, & !< integration point number
el !< element number
2020-07-12 20:14:26 +05:30
real(pReal), intent(in) :: &
2020-12-29 11:50:37 +05:30
T !< plastic velocity gradient
2020-07-12 20:14:26 +05:30
real(pReal), intent(inout) :: &
TDot, &
dTDot_dT
real(pReal) :: &
my_Tdot, &
my_dTdot_dT
integer :: &
2020-12-30 17:15:08 +05:30
ph, &
homog, &
instance, &
2020-12-30 17:15:08 +05:30
me, &
so, &
co
2020-08-13 00:44:00 +05:30
homog = material_homogenizationAt(el)
instance = thermal_typeInstance(homog)
2020-08-13 00:44:00 +05:30
2020-12-30 17:15:08 +05:30
do co = 1, homogenization_Nconstituents(homog)
ph = material_phaseAt(co,el)
me = material_phasememberAt(co,ip,el)
do so = 1, phase_Nsources(ph)
select case(phase_source(so,ph))
case (SOURCE_thermal_dissipation_ID)
call source_thermal_dissipation_getRateAndItsTangent(my_Tdot, my_dTdot_dT, &
2020-12-30 17:15:08 +05:30
mech_S(ph,me),mech_L_p(ph,me), ph)
2020-08-13 00:44:00 +05:30
case (SOURCE_thermal_externalheat_ID)
call source_thermal_externalheat_getRateAndItsTangent(my_Tdot, my_dTdot_dT, &
2020-12-30 17:15:08 +05:30
ph, me)
2020-08-13 00:44:00 +05:30
case default
my_Tdot = 0.0_pReal
my_dTdot_dT = 0.0_pReal
end select
Tdot = Tdot + my_Tdot
dTdot_dT = dTdot_dT + my_dTdot_dT
2020-08-13 00:44:00 +05:30
enddo
enddo
2020-08-13 00:44:00 +05:30
2020-07-12 20:14:26 +05:30
end subroutine constitutive_thermal_getRateAndItsTangents
2020-07-15 18:05:21 +05:30
2020-12-30 14:24:06 +05:30
! getter for non-thermal (e.g. mech)
module function thermal_T(ph,me) result(T)
2020-12-30 14:24:06 +05:30
integer, intent(in) :: ph, me
2020-12-30 14:24:06 +05:30
real(pReal) :: T
T = current(ph)%T(me)
2020-12-30 14:24:06 +05:30
end function thermal_T
2020-12-30 14:24:06 +05:30
2020-12-30 16:30:47 +05:30
! setter for homogenization
module subroutine constitutive_thermal_setT(T,co,ip,el)
real(pReal), intent(in) :: T
integer, intent(in) :: co, ip, el
end subroutine constitutive_thermal_setT
2020-07-15 18:05:21 +05:30
end submodule constitutive_thermal