DAMASK_EICMD/src/thermal_conduction.f90

245 lines
8.7 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @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
use constitutive
2020-08-15 19:32:10 +05:30
use YAML_types
2021-01-08 02:45:18 +05:30
use discretization
2020-03-02 20:19:14 +05:30
implicit none
private
2019-12-21 15:13:36 +05:30
2019-12-21 15:01:19 +05:30
type :: tParameters
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
public :: &
thermal_conduction_init, &
2021-01-11 20:43:59 +05:30
thermal_conduction_getSource, &
2020-02-29 17:27:19 +05:30
thermal_conduction_getConductivity, &
thermal_conduction_getSpecificHeat, &
thermal_conduction_getMassDensity, &
thermal_conduction_putTemperatureAndItsRate, &
2019-12-11 04:40:02 +05:30
thermal_conduction_results
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2021-01-08 02:45:18 +05:30
subroutine thermal_conduction_init(T)
2019-06-07 00:44:37 +05:30
2021-01-08 02:45:18 +05:30
real(pReal), dimension(:), intent(inout) :: T
integer :: Ninstances,Nmaterialpoints,ho,ip,el,ce
2020-08-15 19:32:10 +05:30
class(tNode), pointer :: &
material_homogenization, &
homog, &
homogThermal
2021-01-08 02:45:18 +05:30
print'(/,a)', ' <<<+- thermal_conduction init -+>>>'; flush(6)
2020-03-02 20:19:14 +05:30
Ninstances = count(thermal_type == THERMAL_conduction_ID)
allocate(param(Ninstances))
2020-03-02 20:19:14 +05:30
material_homogenization => config_material%get('homogenization')
2021-01-08 02:45:18 +05:30
do ho = 1, size(material_name_homogenization)
if (thermal_type(ho) /= THERMAL_conduction_ID) cycle
homog => material_homogenization%get(ho)
2020-08-15 19:32:10 +05:30
homogThermal => homog%get('thermal')
2021-01-08 02:45:18 +05:30
associate(prm => param(thermal_typeInstance(ho)))
2020-03-02 20:19:14 +05:30
2020-08-15 19:32:10 +05:30
#if defined (__GFORTRAN__)
prm%output = output_asStrings(homogThermal)
#else
prm%output = homogThermal%get_asStrings('output',defaultVal=emptyStringArray)
#endif
2020-03-02 20:19:14 +05:30
2021-01-08 02:45:18 +05:30
Nmaterialpoints=count(material_homogenizationAt==ho)
2020-03-02 20:19:14 +05:30
2021-01-08 02:45:18 +05:30
allocate (temperature (ho)%p(Nmaterialpoints), source=thermal_initialT(ho))
allocate (temperatureRate(ho)%p(Nmaterialpoints), 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
2021-01-08 02:45:18 +05:30
ce = 0
do el = 1, discretization_Nelems
do ip = 1, discretization_nIPs
ce = ce + 1
ho = material_homogenizationAt(el)
if (thermal_type(ho) == THERMAL_conduction_ID) T(ce) = thermal_initialT(ho)
enddo
enddo
end subroutine thermal_conduction_init
2020-01-31 11:25:26 +05:30
!--------------------------------------------------------------------------------------------------
2020-06-26 15:14:17 +05:30
!> @brief return heat generation rate
!--------------------------------------------------------------------------------------------------
2021-01-11 20:43:59 +05:30
subroutine thermal_conduction_getSource(Tdot, T,ip,el)
2020-03-02 20:19:14 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: &
T
real(pReal), intent(out) :: &
2021-01-08 13:27:30 +05:30
Tdot
integer :: &
homog
2021-01-08 02:45:18 +05:30
2021-01-08 12:56:17 +05:30
homog = material_homogenizationAt(el)
2021-01-11 20:43:59 +05:30
call constitutive_thermal_getRate(TDot, T,ip,el)
2020-03-02 20:19:14 +05:30
Tdot = Tdot/real(homogenization_Nconstituents(homog),pReal)
2020-03-02 20:19:14 +05:30
2021-01-11 20:43:59 +05:30
end subroutine thermal_conduction_getSource
2020-03-02 20:19:14 +05:30
2018-12-31 02:24:50 +05:30
!--------------------------------------------------------------------------------------------------
2020-06-26 15:14:17 +05:30
!> @brief return homogenized thermal conductivity in reference configuration
!--------------------------------------------------------------------------------------------------
2020-02-29 17:27:19 +05:30
function thermal_conduction_getConductivity(ip,el)
2020-03-02 20:19:14 +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
2021-01-08 02:45:18 +05:30
integer :: &
2020-12-29 22:57:24 +05:30
co
2020-03-02 20:19:14 +05:30
2020-02-29 17:27:19 +05:30
thermal_conduction_getConductivity = 0.0_pReal
2021-01-08 02:45:18 +05:30
2020-12-29 22:57:24 +05:30
do co = 1, homogenization_Nconstituents(material_homogenizationAt(el))
2020-02-29 17:27:19 +05:30
thermal_conduction_getConductivity = thermal_conduction_getConductivity + &
2020-12-29 22:57:24 +05:30
crystallite_push33ToRef(co,ip,el,lattice_K(:,:,material_phaseAt(co,el)))
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_Nconstituents(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
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized specific heat capacity
!--------------------------------------------------------------------------------------------------
function thermal_conduction_getSpecificHeat(ip,el)
2020-03-02 20:19:14 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal) :: &
thermal_conduction_getSpecificHeat
2020-12-29 22:57:24 +05:30
integer :: &
2020-12-29 22:57:24 +05:30
co
2020-03-02 20:19:14 +05:30
thermal_conduction_getSpecificHeat = 0.0_pReal
2020-03-02 20:19:14 +05:30
2020-12-29 22:57:24 +05:30
do co = 1, homogenization_Nconstituents(material_homogenizationAt(el))
2020-01-31 11:25:26 +05:30
thermal_conduction_getSpecificHeat = thermal_conduction_getSpecificHeat &
2020-12-29 22:57:24 +05:30
+ lattice_c_p(material_phaseAt(co,el))
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_Nconstituents(material_homogenizationAt(el)),pReal)
2020-03-02 20:19:14 +05:30
end function thermal_conduction_getSpecificHeat
2019-12-21 14:40:23 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized mass density
!--------------------------------------------------------------------------------------------------
function thermal_conduction_getMassDensity(ip,el)
2019-06-07 00:44:37 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal) :: &
thermal_conduction_getMassDensity
2021-01-08 02:45:18 +05:30
integer :: &
2020-12-29 22:57:24 +05:30
co
2020-03-02 20:19:14 +05:30
2020-12-29 22:57:24 +05:30
thermal_conduction_getMassDensity = 0.0_pReal
2020-03-02 20:19:14 +05:30
2020-12-29 22:57:24 +05:30
do co = 1, homogenization_Nconstituents(material_homogenizationAt(el))
thermal_conduction_getMassDensity = thermal_conduction_getMassDensity &
2020-12-29 22:57:24 +05:30
+ lattice_rho(material_phaseAt(co,el))
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_Nconstituents(material_homogenizationAt(el)),pReal)
2020-03-02 20:19:14 +05:30
end function thermal_conduction_getMassDensity
2018-12-31 02:24:50 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief updates thermal state with solution from heat conduction PDE
!--------------------------------------------------------------------------------------------------
subroutine thermal_conduction_putTemperatureAndItsRate(T,Tdot,ip,el)
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
homog = material_homogenizationAt(el)
offset = material_homogenizationMemberAt(ip,el)
temperature (homog)%p(offset) = T
temperatureRate(homog)%p(offset) = Tdot
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
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)))
outputsLoop: do o = 1,size(prm%output)
select case(trim(prm%output(o)))
2020-08-15 19:32:10 +05:30
case('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
end module thermal_conduction