2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @brief material subroutine for temperature evolution from heat conduction
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
module thermal_conduction
|
2019-06-07 00:44:37 +05:30
|
|
|
use prec
|
|
|
|
use material
|
|
|
|
use config
|
|
|
|
use lattice
|
2019-12-11 00:55:19 +05:30
|
|
|
use results
|
2019-06-07 00:44:37 +05:30
|
|
|
use crystallite
|
2020-07-09 04:31:08 +05:30
|
|
|
use constitutive
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
implicit none
|
|
|
|
private
|
2019-12-21 15:13:36 +05:30
|
|
|
|
2019-12-21 15:01:19 +05:30
|
|
|
type :: tParameters
|
2020-02-15 03:20:30 +05:30
|
|
|
character(len=pStringLen), allocatable, dimension(:) :: &
|
|
|
|
output
|
2019-12-21 15:01:19 +05:30
|
|
|
end type tParameters
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-12-21 15:01:19 +05:30
|
|
|
type(tparameters), dimension(:), allocatable :: &
|
2020-03-02 20:19:14 +05:30
|
|
|
param
|
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
public :: &
|
|
|
|
thermal_conduction_init, &
|
|
|
|
thermal_conduction_getSourceAndItsTangent, &
|
2020-02-29 17:27:19 +05:30
|
|
|
thermal_conduction_getConductivity, &
|
2019-03-24 16:37:51 +05:30
|
|
|
thermal_conduction_getSpecificHeat, &
|
|
|
|
thermal_conduction_getMassDensity, &
|
|
|
|
thermal_conduction_putTemperatureAndItsRate, &
|
2019-12-11 04:40:02 +05:30
|
|
|
thermal_conduction_results
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief module initialization
|
|
|
|
!> @details reads in material parameters, allocates arrays, and does sanity checks
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-12-31 03:00:21 +05:30
|
|
|
subroutine thermal_conduction_init
|
2019-06-07 00:44:37 +05:30
|
|
|
|
2020-03-02 20:19:14 +05:30
|
|
|
integer :: Ninstance,NofMyHomog,h
|
|
|
|
|
2020-02-15 03:20:30 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- thermal_'//THERMAL_CONDUCTION_label//' init -+>>>'; flush(6)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
|
|
|
Ninstance = count(thermal_type == THERMAL_conduction_ID)
|
|
|
|
allocate(param(Ninstance))
|
|
|
|
|
|
|
|
do h = 1, size(config_homogenization)
|
2019-12-21 15:25:11 +05:30
|
|
|
if (thermal_type(h) /= THERMAL_conduction_ID) cycle
|
|
|
|
associate(prm => param(thermal_typeInstance(h)),config => config_homogenization(h))
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-02-15 03:20:30 +05:30
|
|
|
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-12-21 15:25:11 +05:30
|
|
|
NofMyHomog=count(material_homogenizationAt==h)
|
|
|
|
thermalState(h)%sizeState = 0
|
|
|
|
allocate(thermalState(h)%state0 (0,NofMyHomog))
|
|
|
|
allocate(thermalState(h)%subState0(0,NofMyHomog))
|
|
|
|
allocate(thermalState(h)%state (0,NofMyHomog))
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-01-23 17:18:18 +05:30
|
|
|
thermalMapping(h)%p => material_homogenizationMemberAt
|
2019-12-21 15:25:11 +05:30
|
|
|
deallocate(temperature (h)%p)
|
|
|
|
allocate (temperature (h)%p(NofMyHomog), source=thermal_initialT(h))
|
|
|
|
deallocate(temperatureRate(h)%p)
|
|
|
|
allocate (temperatureRate(h)%p(NofMyHomog), source=0.0_pReal)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-12-21 15:01:19 +05:30
|
|
|
end associate
|
2019-12-21 15:25:11 +05:30
|
|
|
enddo
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end subroutine thermal_conduction_init
|
|
|
|
|
2020-01-31 11:25:26 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-06-26 15:14:17 +05:30
|
|
|
!> @brief return heat generation rate
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine thermal_conduction_getSourceAndItsTangent(Tdot, dTdot_dT, T, ip, el)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
T
|
|
|
|
real(pReal), intent(out) :: &
|
|
|
|
Tdot, dTdot_dT
|
|
|
|
integer :: &
|
2020-07-09 04:31:08 +05:30
|
|
|
homog
|
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
Tdot = 0.0_pReal
|
|
|
|
dTdot_dT = 0.0_pReal
|
2020-07-09 04:31:08 +05:30
|
|
|
|
|
|
|
homog = material_homogenizationAt(el)
|
|
|
|
call constitutive_thermal_getRateAndItsTangents(TDot, dTDot_dT, T, crystallite_S,crystallite_Lp ,ip, el)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
Tdot = Tdot/real(homogenization_Ngrains(homog),pReal)
|
|
|
|
dTdot_dT = dTdot_dT/real(homogenization_Ngrains(homog),pReal)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end subroutine thermal_conduction_getSourceAndItsTangent
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2018-12-31 02:24:50 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-06-26 15:14:17 +05:30
|
|
|
!> @brief return homogenized thermal conductivity in reference configuration
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-29 17:27:19 +05:30
|
|
|
function thermal_conduction_getConductivity(ip,el)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), dimension(3,3) :: &
|
2020-02-29 17:27:19 +05:30
|
|
|
thermal_conduction_getConductivity
|
2019-03-24 16:37:51 +05:30
|
|
|
integer :: &
|
|
|
|
grain
|
2020-03-02 20:19:14 +05:30
|
|
|
|
|
|
|
|
2020-02-29 17:27:19 +05:30
|
|
|
thermal_conduction_getConductivity = 0.0_pReal
|
2019-06-07 00:44:37 +05:30
|
|
|
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
|
2020-02-29 17:27:19 +05:30
|
|
|
thermal_conduction_getConductivity = thermal_conduction_getConductivity + &
|
|
|
|
crystallite_push33ToRef(grain,ip,el,lattice_thermalConductivity(:,:,material_phaseAt(grain,el)))
|
2019-03-24 16:37:51 +05:30
|
|
|
enddo
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-02-29 17:27:19 +05:30
|
|
|
thermal_conduction_getConductivity = thermal_conduction_getConductivity &
|
|
|
|
/ real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-02-29 17:27:19 +05:30
|
|
|
end function thermal_conduction_getConductivity
|
2018-12-31 02:24:50 +05:30
|
|
|
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns homogenized specific heat capacity
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function thermal_conduction_getSpecificHeat(ip,el)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal) :: &
|
|
|
|
thermal_conduction_getSpecificHeat
|
|
|
|
integer :: &
|
|
|
|
grain
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
thermal_conduction_getSpecificHeat = 0.0_pReal
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-06-07 00:44:37 +05:30
|
|
|
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
|
2020-01-31 11:25:26 +05:30
|
|
|
thermal_conduction_getSpecificHeat = thermal_conduction_getSpecificHeat &
|
|
|
|
+ lattice_specificHeat(material_phaseAt(grain,el))
|
2019-03-24 16:37:51 +05:30
|
|
|
enddo
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-01-31 11:25:26 +05:30
|
|
|
thermal_conduction_getSpecificHeat = thermal_conduction_getSpecificHeat &
|
|
|
|
/ real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end function thermal_conduction_getSpecificHeat
|
2019-12-21 14:40:23 +05:30
|
|
|
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief returns homogenized mass density
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
function thermal_conduction_getMassDensity(ip,el)
|
2019-06-07 00:44:37 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal) :: &
|
|
|
|
thermal_conduction_getMassDensity
|
|
|
|
integer :: &
|
|
|
|
grain
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
thermal_conduction_getMassDensity = 0.0_pReal
|
2020-03-02 20:19:14 +05:30
|
|
|
|
|
|
|
|
2019-06-07 00:44:37 +05:30
|
|
|
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
|
2019-03-24 16:37:51 +05:30
|
|
|
thermal_conduction_getMassDensity = thermal_conduction_getMassDensity &
|
2019-06-15 17:27:24 +05:30
|
|
|
+ lattice_massDensity(material_phaseAt(grain,el))
|
2019-03-24 16:37:51 +05:30
|
|
|
enddo
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2020-01-31 11:25:26 +05:30
|
|
|
thermal_conduction_getMassDensity = thermal_conduction_getMassDensity &
|
|
|
|
/ real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end function thermal_conduction_getMassDensity
|
2018-12-31 02:24:50 +05:30
|
|
|
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief updates thermal state with solution from heat conduction PDE
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine thermal_conduction_putTemperatureAndItsRate(T,Tdot,ip,el)
|
2019-03-24 16:37:51 +05:30
|
|
|
|
|
|
|
integer, intent(in) :: &
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(in) :: &
|
|
|
|
T, &
|
|
|
|
Tdot
|
|
|
|
integer :: &
|
|
|
|
homog, &
|
2020-03-02 20:19:14 +05:30
|
|
|
offset
|
|
|
|
|
2019-03-24 16:37:51 +05:30
|
|
|
homog = material_homogenizationAt(el)
|
|
|
|
offset = thermalMapping(homog)%p(ip,el)
|
|
|
|
temperature (homog)%p(offset) = T
|
|
|
|
temperatureRate(homog)%p(offset) = Tdot
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
end subroutine thermal_conduction_putTemperatureAndItsRate
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-12-11 00:55:19 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief writes results to HDF5 output file
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine thermal_conduction_results(homog,group)
|
|
|
|
|
|
|
|
integer, intent(in) :: homog
|
|
|
|
character(len=*), intent(in) :: group
|
2020-02-15 03:20:30 +05:30
|
|
|
|
2019-12-21 15:13:36 +05:30
|
|
|
integer :: o
|
2020-03-02 20:19:14 +05:30
|
|
|
|
2019-12-21 15:13:36 +05:30
|
|
|
associate(prm => param(damage_typeInstance(homog)))
|
2020-02-15 03:20:30 +05:30
|
|
|
outputsLoop: do o = 1,size(prm%output)
|
|
|
|
select case(trim(prm%output(o)))
|
|
|
|
case('temperature') ! ToDo: should be 'T'
|
2019-12-11 00:55:19 +05:30
|
|
|
call results_writeDataset(group,temperature(homog)%p,'T',&
|
|
|
|
'temperature','K')
|
|
|
|
end select
|
|
|
|
enddo outputsLoop
|
2019-12-21 15:13:36 +05:30
|
|
|
end associate
|
2019-12-11 00:55:19 +05:30
|
|
|
|
|
|
|
end subroutine thermal_conduction_results
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end module thermal_conduction
|