DAMASK_EICMD/src/thermal_adiabatic.f90

268 lines
9.8 KiB
Fortran
Raw Normal View History

2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine for adiabatic temperature evolution
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
module thermal_adiabatic
use prec
use config
use numerics
use material
2019-12-11 00:55:19 +05:30
use results
use source_thermal_dissipation
use source_thermal_externalheat
use crystallite
use lattice
2019-03-24 15:58:20 +05:30
implicit none
private
2019-12-21 15:13:36 +05:30
2019-03-24 15:58:20 +05:30
enum, bind(c)
enumerator :: undefined_ID, &
temperature_ID
end enum
2019-12-21 15:01:19 +05:30
type :: tParameters
integer(kind(undefined_ID)), dimension(:), allocatable :: &
outputID
end type tParameters
type(tparameters), dimension(:), allocatable :: &
param
2019-03-24 15:58:20 +05:30
public :: &
thermal_adiabatic_init, &
thermal_adiabatic_updateState, &
thermal_adiabatic_getSourceAndItsTangent, &
thermal_adiabatic_getSpecificHeat, &
thermal_adiabatic_getMassDensity, &
2019-12-11 04:40:02 +05:30
thermal_adiabatic_results
2019-03-24 15:58:20 +05:30
2014-10-09 19:36:45 +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_adiabatic_init
2014-10-09 19:36:45 +05:30
2019-12-21 15:25:11 +05:30
integer :: maxNinstance,o,h,NofMyHomog
character(len=pStringLen), dimension(:), allocatable :: outputs
2014-10-09 19:36:45 +05:30
2019-12-21 15:01:19 +05:30
write(6,'(/,a)') ' <<<+- thermal_'//THERMAL_ADIABATIC_label//' init -+>>>'; flush(6)
2019-03-24 15:58:20 +05:30
maxNinstance = count(thermal_type == THERMAL_adiabatic_ID)
if (maxNinstance == 0) return
2019-12-21 15:01:19 +05:30
allocate(param(maxNinstance))
2019-12-21 15:25:11 +05:30
do h = 1, size(thermal_type)
if (thermal_type(h) /= THERMAL_adiabatic_ID) cycle
associate(prm => param(thermal_typeInstance(h)),config => config_homogenization(h))
2019-12-21 15:01:19 +05:30
2019-12-21 15:13:36 +05:30
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
allocate(prm%outputID(0))
2019-12-21 15:25:11 +05:30
do o=1, size(outputs)
select case(outputs(o))
2019-03-24 15:58:20 +05:30
case('temperature')
2019-12-21 15:25:11 +05:30
prm%outputID = [prm%outputID, temperature_ID]
2019-03-24 15:58:20 +05:30
end select
enddo
2019-12-21 15:25:11 +05:30
NofMyHomog=count(material_homogenizationAt==h)
thermalState(h)%sizeState = 1
allocate(thermalState(h)%state0 (1,NofMyHomog), source=thermal_initialT(h))
allocate(thermalState(h)%subState0(1,NofMyHomog), source=thermal_initialT(h))
allocate(thermalState(h)%state (1,NofMyHomog), source=thermal_initialT(h))
2019-03-24 15:58:20 +05:30
2019-12-21 15:25:11 +05:30
nullify(thermalMapping(h)%p)
thermalMapping(h)%p => material_homogenizationMemberAt
2019-12-21 15:25:11 +05:30
deallocate(temperature(h)%p)
temperature(h)%p => thermalState(h)%state(1,:)
deallocate(temperatureRate(h)%p)
allocate (temperatureRate(h)%p(NofMyHomog), source=0.0_pReal)
2019-12-21 15:01:19 +05:30
end associate
2019-12-21 15:25:11 +05:30
enddo
2018-12-31 03:00:21 +05:30
2014-10-09 19:36:45 +05:30
end subroutine thermal_adiabatic_init
2019-03-24 15:58:20 +05:30
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates adiabatic change in temperature based on local heat generation model
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
function thermal_adiabatic_updateState(subdt, ip, el)
2019-03-24 15:58:20 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: &
subdt
logical, dimension(2) :: &
thermal_adiabatic_updateState
integer :: &
homog, &
offset
real(pReal) :: &
T, Tdot, dTdot_dT
homog = material_homogenizationAt(el)
offset = material_homogenizationMemberAt(ip,el)
2019-03-24 15:58:20 +05:30
T = thermalState(homog)%subState0(1,offset)
call thermal_adiabatic_getSourceAndItsTangent(Tdot, dTdot_dT, T, ip, el)
T = T + subdt*Tdot/(thermal_adiabatic_getSpecificHeat(ip,el)*thermal_adiabatic_getMassDensity(ip,el))
thermal_adiabatic_updateState = [ abs(T - thermalState(homog)%state(1,offset)) &
<= err_thermal_tolAbs &
.or. abs(T - thermalState(homog)%state(1,offset)) &
<= err_thermal_tolRel*abs(thermalState(homog)%state(1,offset)), &
.true.]
2019-03-24 15:58:20 +05:30
temperature (homog)%p(thermalMapping(homog)%p(ip,el)) = T
temperatureRate(homog)%p(thermalMapping(homog)%p(ip,el)) = &
(thermalState(homog)%state(1,offset) - thermalState(homog)%subState0(1,offset))/(subdt+tiny(0.0_pReal))
end function thermal_adiabatic_updateState
!--------------------------------------------------------------------------------------------------
!> @brief returns heat generation rate
!--------------------------------------------------------------------------------------------------
subroutine thermal_adiabatic_getSourceAndItsTangent(Tdot, dTdot_dT, T, ip, el)
2019-03-24 15:58:20 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), intent(in) :: &
T
real(pReal), intent(out) :: &
Tdot, dTdot_dT
real(pReal) :: &
my_Tdot, my_dTdot_dT
integer :: &
phase, &
homog, &
instance, &
grain, &
source, &
constituent
homog = material_homogenizationAt(el)
instance = thermal_typeInstance(homog)
2019-03-24 15:58:20 +05:30
Tdot = 0.0_pReal
dTdot_dT = 0.0_pReal
do grain = 1, homogenization_Ngrains(homog)
phase = material_phaseAt(grain,el)
constituent = material_phasememberAt(grain,ip,el)
2019-03-24 15:58:20 +05:30
do source = 1, phase_Nsources(phase)
select case(phase_source(source,phase))
case (SOURCE_thermal_dissipation_ID)
call source_thermal_dissipation_getRateAndItsTangent(my_Tdot, my_dTdot_dT, &
crystallite_S(1:3,1:3,grain,ip,el), &
crystallite_Lp(1:3,1:3,grain,ip,el), &
phase)
2014-10-09 19:36:45 +05:30
2019-03-24 15:58:20 +05:30
case (SOURCE_thermal_externalheat_ID)
call source_thermal_externalheat_getRateAndItsTangent(my_Tdot, my_dTdot_dT, &
phase, constituent)
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
enddo
enddo
Tdot = Tdot/real(homogenization_Ngrains(homog),pReal)
dTdot_dT = dTdot_dT/real(homogenization_Ngrains(homog),pReal)
2014-10-09 19:36:45 +05:30
end subroutine thermal_adiabatic_getSourceAndItsTangent
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized specific heat capacity
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
function thermal_adiabatic_getSpecificHeat(ip,el)
2019-03-24 15:58:20 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
2019-03-24 15:58:20 +05:30
real(pReal) :: &
thermal_adiabatic_getSpecificHeat
integer :: &
grain
thermal_adiabatic_getSpecificHeat = 0.0_pReal
2019-06-07 00:44:37 +05:30
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
2019-03-24 15:58:20 +05:30
thermal_adiabatic_getSpecificHeat = thermal_adiabatic_getSpecificHeat + &
2019-06-15 17:27:24 +05:30
lattice_specificHeat(material_phaseAt(grain,el))
2019-03-24 15:58:20 +05:30
enddo
2019-03-24 15:58:20 +05:30
thermal_adiabatic_getSpecificHeat = &
2019-06-07 00:44:37 +05:30
thermal_adiabatic_getSpecificHeat/real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
2019-03-24 15:58:20 +05:30
end function thermal_adiabatic_getSpecificHeat
2018-12-31 02:24:50 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized mass density
!--------------------------------------------------------------------------------------------------
function thermal_adiabatic_getMassDensity(ip,el)
2019-03-24 15:58:20 +05:30
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal) :: &
thermal_adiabatic_getMassDensity
integer :: &
grain
2019-03-24 15:58:20 +05:30
thermal_adiabatic_getMassDensity = 0.0_pReal
2019-06-07 00:44:37 +05:30
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
2019-03-24 15:58:20 +05:30
thermal_adiabatic_getMassDensity = thermal_adiabatic_getMassDensity + &
2019-06-15 17:27:24 +05:30
lattice_massDensity(material_phaseAt(grain,el))
2019-03-24 15:58:20 +05:30
enddo
thermal_adiabatic_getMassDensity = &
2019-06-07 00:44:37 +05:30
thermal_adiabatic_getMassDensity/real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
end function thermal_adiabatic_getMassDensity
2018-12-31 02:24:50 +05:30
2019-12-11 00:55:19 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine thermal_adiabatic_results(homog,group)
integer, intent(in) :: homog
character(len=*), intent(in) :: group
2019-12-21 15:13:36 +05:30
integer :: o
2019-12-11 00:55:19 +05:30
2019-12-21 15:13:36 +05:30
associate(prm => param(damage_typeInstance(homog)))
2019-12-11 00:55:19 +05:30
2019-12-21 15:13:36 +05:30
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
2019-12-11 00:55:19 +05:30
case (temperature_ID)
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_adiabatic_results
2014-10-09 19:36:45 +05:30
end module thermal_adiabatic