DAMASK_EICMD/src/thermal_adiabatic.f90

230 lines
8.6 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 material
2019-12-11 00:55:19 +05:30
use results
use constitutive
2020-08-15 19:32:10 +05:30
use YAML_types
use crystallite
use lattice
2019-03-24 15:58:20 +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
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
2020-08-15 19:32:10 +05:30
integer :: maxNinstance,h,NofMyHomog
class(tNode), pointer :: &
material_homogenization, &
homog, &
homogThermal
print'(/,a)', ' <<<+- thermal_adiabatic 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))
material_homogenization => config_material%get('homogenization')
2020-08-15 19:32:10 +05:30
do h = 1, material_Nhomogenization
2019-12-21 15:25:11 +05:30
if (thermal_type(h) /= THERMAL_adiabatic_ID) cycle
2020-08-15 19:32:10 +05:30
homog => material_homogenization%get(h)
homogThermal => homog%get('thermal')
associate(prm => param(thermal_typeInstance(h)))
#if defined (__GFORTRAN__)
prm%output = output_asStrings(homogThermal)
#else
prm%output = homogThermal%get_asStrings('output',defaultVal=emptyStringArray)
#endif
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
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)) &
2020-06-16 21:23:14 +05:30
<= 1.0e-2_pReal &
2019-03-24 15:58:20 +05:30
.or. abs(T - thermalState(homog)%state(1,offset)) &
2020-06-16 21:23:14 +05:30
<= 1.0e-6_pReal*abs(thermalState(homog)%state(1,offset)), &
2019-03-24 15:58:20 +05:30
.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
integer :: &
homog
2019-03-24 15:58:20 +05:30
Tdot = 0.0_pReal
dTdot_dT = 0.0_pReal
2014-10-09 19:36:45 +05:30
homog = material_homogenizationAt(el)
2020-07-15 18:05:21 +05:30
call constitutive_thermal_getRateAndItsTangents(TDot, dTDot_dT, T, crystallite_S, crystallite_Lp, ip, el)
2019-03-24 15:58:20 +05:30
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))
2020-01-31 11:22:11 +05:30
thermal_adiabatic_getSpecificHeat = thermal_adiabatic_getSpecificHeat &
+ lattice_specificHeat(material_phaseAt(grain,el))
2019-03-24 15:58:20 +05:30
enddo
2020-01-31 11:22:11 +05:30
thermal_adiabatic_getSpecificHeat = 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))
2020-01-31 11:22:11 +05:30
thermal_adiabatic_getMassDensity = thermal_adiabatic_getMassDensity &
+ lattice_massDensity(material_phaseAt(grain,el))
2019-03-24 15:58:20 +05:30
enddo
2020-01-31 11:22:11 +05:30
thermal_adiabatic_getMassDensity = 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)))
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_adiabatic_results
2014-10-09 19:36:45 +05:30
end module thermal_adiabatic