DAMASK_EICMD/src/phase_thermal_externalheat.f90

138 lines
5.5 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Michigan State University
!> @brief material subroutine for variable heat source
!--------------------------------------------------------------------------------------------------
submodule(phase:thermal) externalheat
integer, dimension(:), allocatable :: &
2019-12-05 15:50:05 +05:30
source_thermal_externalheat_offset, & !< which source is my current thermal dissipation mechanism?
source_thermal_externalheat_instance !< instance of thermal dissipation source mechanism
2020-08-15 19:32:10 +05:30
type :: tParameters !< container type for internal constitutive parameters
real(pReal), dimension(:), allocatable :: &
2021-01-08 02:45:18 +05:30
t_n, &
2020-09-23 04:38:13 +05:30
f_T
integer :: &
nIntervals
end type tParameters
type(tParameters), dimension(:), allocatable :: param !< containers of constitutive parameters (len Ninstances)
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2021-01-26 12:25:06 +05:30
module function externalheat_init(source_length) result(mySources)
2020-02-26 23:07:17 +05:30
2021-01-08 02:45:18 +05:30
integer, intent(in) :: source_length
2020-08-15 19:32:10 +05:30
logical, dimension(:,:), allocatable :: mySources
class(tNode), pointer :: &
phases, &
phase, &
2021-01-08 02:45:18 +05:30
sources, thermal, &
src
integer :: Ninstances,sourceOffset,Nconstituents,p
2020-02-26 23:07:17 +05:30
2021-01-08 12:07:51 +05:30
print'(/,a)', ' <<<+- thermal_externalheat init -+>>>'
2020-02-26 23:07:17 +05:30
2021-01-08 02:45:18 +05:30
mySources = thermal_active('externalheat',source_length)
Ninstances = count(mySources)
print'(a,i2)', ' # instances: ',Ninstances; flush(IO_STDOUT)
if(Ninstances == 0) return
2020-02-26 23:07:17 +05:30
phases => config_material%get('phase')
allocate(param(Ninstances))
2020-08-15 19:32:10 +05:30
allocate(source_thermal_externalheat_offset (phases%length), source=0)
allocate(source_thermal_externalheat_instance(phases%length), source=0)
do p = 1, phases%length
2021-01-08 02:45:18 +05:30
phase => phases%get(p)
2020-08-15 19:32:10 +05:30
if(any(mySources(:,p))) source_thermal_externalheat_instance(p) = count(mySources(:,1:p))
if(count(mySources(:,p)) == 0) cycle
2021-01-08 02:45:18 +05:30
thermal => phase%get('thermal')
sources => thermal%get('source')
2020-08-15 19:32:10 +05:30
do sourceOffset = 1, sources%length
if(mySources(sourceOffset,p)) then
2020-02-29 12:33:06 +05:30
source_thermal_externalheat_offset(p) = sourceOffset
2020-08-15 19:32:10 +05:30
associate(prm => param(source_thermal_externalheat_instance(p)))
2021-01-08 02:45:18 +05:30
src => sources%get(sourceOffset)
2020-02-26 23:07:17 +05:30
2020-09-23 04:38:13 +05:30
prm%t_n = src%get_asFloats('t_n')
prm%nIntervals = size(prm%t_n) - 1
2020-02-26 23:07:17 +05:30
2020-09-23 04:38:13 +05:30
prm%f_T = src%get_asFloats('f_T',requiredSize = size(prm%t_n))
2020-02-26 23:07:17 +05:30
Nconstituents = count(material_phaseAt==p) * discretization_nIPs
2021-01-08 02:45:18 +05:30
call constitutive_allocateState(thermalState(p)%p(sourceOffset),Nconstituents,1,1,0)
2020-08-15 19:32:10 +05:30
end associate
endif
enddo
enddo
2021-01-26 12:25:06 +05:30
end function externalheat_init
2019-02-22 19:51:48 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief rate of change of state
!> @details state only contains current time to linearly interpolate given heat powers
!--------------------------------------------------------------------------------------------------
2021-01-26 12:25:06 +05:30
module subroutine externalheat_dotState(ph, me)
2020-02-26 23:07:17 +05:30
integer, intent(in) :: &
2021-01-19 15:02:56 +05:30
ph, &
me
2020-02-29 12:28:33 +05:30
integer :: &
sourceOffset
2020-02-26 23:07:17 +05:30
2021-01-19 15:02:56 +05:30
sourceOffset = source_thermal_externalheat_offset(ph)
2020-02-26 23:07:17 +05:30
2021-01-19 15:02:56 +05:30
thermalState(ph)%p(sourceOffset)%dotState(1,me) = 1.0_pReal ! state is current time
2021-01-26 12:25:06 +05:30
end subroutine externalheat_dotState
2019-12-05 15:50:05 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-26 23:07:17 +05:30
!> @brief returns local heat generation rate
!--------------------------------------------------------------------------------------------------
2021-01-19 15:02:56 +05:30
module subroutine thermal_externalheat_getRate(TDot, ph, me)
2020-02-26 23:07:17 +05:30
integer, intent(in) :: &
2021-01-19 15:02:56 +05:30
ph, &
me
real(pReal), intent(out) :: &
2021-01-08 13:27:30 +05:30
TDot
2020-02-29 12:28:33 +05:30
integer :: &
2020-02-29 12:28:33 +05:30
sourceOffset, interval
real(pReal) :: &
2020-02-26 23:07:17 +05:30
frac_time
2021-01-19 15:02:56 +05:30
sourceOffset = source_thermal_externalheat_offset(ph)
2020-02-26 23:07:17 +05:30
2021-01-19 15:02:56 +05:30
associate(prm => param(source_thermal_externalheat_instance(ph)))
2020-02-29 19:04:19 +05:30
do interval = 1, prm%nIntervals ! scan through all rate segments
2021-01-19 15:02:56 +05:30
frac_time = (thermalState(ph)%p(sourceOffset)%state(1,me) - prm%t_n(interval)) &
2020-09-23 04:38:13 +05:30
/ (prm%t_n(interval+1) - prm%t_n(interval)) ! fractional time within segment
if ( (frac_time < 0.0_pReal .and. interval == 1) &
2020-02-29 12:28:33 +05:30
.or. (frac_time >= 1.0_pReal .and. interval == prm%nIntervals) &
.or. (frac_time >= 0.0_pReal .and. frac_time < 1.0_pReal) ) &
2020-09-23 04:38:13 +05:30
TDot = prm%f_T(interval ) * (1.0_pReal - frac_time) + &
prm%f_T(interval+1) * frac_time ! interpolate heat rate between segment boundaries...
2021-01-19 15:02:56 +05:30
! ...or extrapolate if outside me bounds
enddo
2020-02-29 12:28:33 +05:30
end associate
2020-02-26 23:07:17 +05:30
2021-01-11 20:43:59 +05:30
end subroutine thermal_externalheat_getRate
2020-02-26 23:07:17 +05:30
2021-01-26 05:50:45 +05:30
end submodule externalheat