helper routine to get heat generation rate needed for MARC/Abaqus
This commit is contained in:
parent
d51a28a2d9
commit
3f14ebe43d
|
@ -43,6 +43,7 @@ module constitutive
|
|||
constitutive_getAdiabaticTemperature, &
|
||||
constitutive_putAdiabaticTemperature, &
|
||||
constitutive_getTemperature, &
|
||||
constitutive_getHeatGeneration, &
|
||||
constitutive_getLocalVacancyConcentration, &
|
||||
constitutive_putLocalVacancyConcentration, &
|
||||
constitutive_getVacancyConcentration, &
|
||||
|
@ -1607,6 +1608,41 @@ function constitutive_getTemperature(ipc, ip, el)
|
|||
|
||||
end function constitutive_getTemperature
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns heat generation rate
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
function constitutive_getHeatGeneration(Tstar_v, Lp, ipc, ip, el)
|
||||
use prec, only: &
|
||||
pReal
|
||||
use material, only: &
|
||||
material_phase, &
|
||||
LOCAL_THERMAL_isothermal_ID, &
|
||||
LOCAL_THERMAL_adiabatic_ID, &
|
||||
phase_thermal
|
||||
use thermal_adiabatic, only: &
|
||||
thermal_adiabatic_getHeatGeneration
|
||||
|
||||
implicit none
|
||||
integer(pInt), intent(in) :: &
|
||||
ipc, & !< grain number
|
||||
ip, & !< integration point number
|
||||
el !< element number
|
||||
real(pReal), intent(in), dimension(6) :: &
|
||||
Tstar_v !< 2nd Piola-Kirchhoff stress
|
||||
real(pReal), intent(in), dimension(3,3) :: &
|
||||
Lp !< plastic velocity gradient
|
||||
real(pReal) :: constitutive_getHeatGeneration
|
||||
|
||||
select case (phase_thermal(material_phase(ipc,ip,el)))
|
||||
case (LOCAL_THERMAL_isothermal_ID)
|
||||
constitutive_getHeatGeneration = 0.0_pReal
|
||||
|
||||
case (LOCAL_THERMAL_adiabatic_ID)
|
||||
constitutive_getHeatGeneration = thermal_adiabatic_getHeatGeneration(Tstar_v, Lp)
|
||||
end select
|
||||
|
||||
end function constitutive_getHeatGeneration
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns local vacancy concentration
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -67,6 +67,7 @@ module homogenization
|
|||
field_putFieldDamage, &
|
||||
field_getLocalTemperature, &
|
||||
field_putFieldTemperature, &
|
||||
field_getHeatGeneration, &
|
||||
field_getLocalVacancyConcentration, &
|
||||
field_putFieldVacancyConcentration, &
|
||||
field_getDamageMobility, &
|
||||
|
@ -1338,6 +1339,38 @@ subroutine field_putFieldTemperature(ip,el,fieldThermalValue)
|
|||
|
||||
end subroutine field_putFieldTemperature
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief return heat generation rate
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
real(pReal) function field_getHeatGeneration(ip,el)
|
||||
use mesh, only: &
|
||||
mesh_element
|
||||
use material, only: &
|
||||
homogenization_Ngrains
|
||||
use crystallite, only: &
|
||||
crystallite_Tstar_v, &
|
||||
crystallite_Lp
|
||||
use constitutive, only: &
|
||||
constitutive_getHeatGeneration
|
||||
|
||||
implicit none
|
||||
integer(pInt), intent(in) :: &
|
||||
ip, & !< integration point number
|
||||
el !< element number
|
||||
integer(pInt) :: &
|
||||
ipc
|
||||
|
||||
field_getHeatGeneration = 0.0_pReal
|
||||
do ipc = 1, homogenization_Ngrains(mesh_element(3,el))
|
||||
field_getHeatGeneration = field_getHeatGeneration + &
|
||||
constitutive_getHeatGeneration(crystallite_Tstar_v(1:6,ipc,ip,el), &
|
||||
crystallite_Lp (1:3,1:3,ipc,ip,el), &
|
||||
ipc,ip,el)
|
||||
enddo
|
||||
field_getHeatGeneration = field_getHeatGeneration/homogenization_Ngrains(mesh_element(3,el))
|
||||
|
||||
end function field_getHeatGeneration
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief ToDo
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -46,6 +46,7 @@ module thermal_adiabatic
|
|||
thermal_adiabatic_getPartionedFT0, &
|
||||
thermal_adiabatic_getTemperature, &
|
||||
thermal_adiabatic_putTemperature, &
|
||||
thermal_adiabatic_getHeatGeneration, &
|
||||
thermal_adiabatic_postResults
|
||||
|
||||
contains
|
||||
|
@ -486,6 +487,25 @@ subroutine thermal_adiabatic_putTemperature(ipc, ip, el, localTemperature)
|
|||
|
||||
end subroutine thermal_adiabatic_putTemperature
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns heat generation rate
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
function thermal_adiabatic_getHeatGeneration(Tstar_v, Lp)
|
||||
use math, only: &
|
||||
math_Mandel6to33
|
||||
|
||||
implicit none
|
||||
real(pReal), intent(in), dimension(6) :: &
|
||||
Tstar_v !< 2nd Piola-Kirchhoff stress
|
||||
real(pReal), intent(in), dimension(3,3) :: &
|
||||
Lp !< plastic velocity gradient
|
||||
real(pReal) :: thermal_adiabatic_getHeatGeneration
|
||||
|
||||
thermal_adiabatic_getHeatGeneration = 0.95_pReal &
|
||||
* sum(abs(math_Mandel6to33(Tstar_v))*Lp)
|
||||
|
||||
end function thermal_adiabatic_getHeatGeneration
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief return array of constitutive results
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
Loading…
Reference in New Issue