2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-02-23 01:07:41 +05:30
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
2015-07-27 16:39:37 +05:30
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
2016-10-19 01:53:52 +05:30
|
|
|
!> @author Philip Eisenlohr, Michigan State University
|
|
|
|
!> @brief material subroutine for variable heat source
|
2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-01-27 01:22:48 +05:30
|
|
|
submodule(phase:thermal) externalheat
|
2015-07-27 16:39:37 +05:30
|
|
|
|
2019-05-28 15:36:21 +05:30
|
|
|
|
2019-12-21 12:25:42 +05:30
|
|
|
integer, dimension(:), allocatable :: &
|
2021-02-13 23:11:30 +05:30
|
|
|
source_thermal_externalheat_offset !< which source is my current thermal dissipation mechanism?
|
2018-12-31 01:28:38 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
type :: tParameters !< container type for internal constitutive parameters
|
2019-02-23 01:07:41 +05:30
|
|
|
real(pReal), dimension(:), allocatable :: &
|
2021-01-08 02:45:18 +05:30
|
|
|
t_n, &
|
2020-09-23 04:38:13 +05:30
|
|
|
f_T
|
2019-03-24 15:48:59 +05:30
|
|
|
integer :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
nIntervals
|
|
|
|
end type tParameters
|
2018-12-31 01:28:38 +05:30
|
|
|
|
2020-10-28 02:03:30 +05:30
|
|
|
type(tParameters), dimension(:), allocatable :: param !< containers of constitutive parameters (len Ninstances)
|
2018-12-31 01:28:38 +05:30
|
|
|
|
|
|
|
|
2015-07-27 16:39:37 +05:30
|
|
|
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
|
2021-02-13 23:11:30 +05:30
|
|
|
integer :: so,Nconstituents,ph
|
2020-02-26 23:07:17 +05:30
|
|
|
|
|
|
|
|
2021-01-08 02:45:18 +05:30
|
|
|
mySources = thermal_active('externalheat',source_length)
|
2021-02-13 23:11:30 +05:30
|
|
|
if(count(mySources) == 0) return
|
|
|
|
print'(/,a)', ' <<<+- phase:thermal:externalheat init -+>>>'
|
|
|
|
print'(a,i2)', ' # phases: ',count(mySources); flush(IO_STDOUT)
|
2021-01-08 02:45:18 +05:30
|
|
|
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2020-09-13 14:09:17 +05:30
|
|
|
phases => config_material%get('phase')
|
2021-01-27 04:26:20 +05:30
|
|
|
allocate(param(phases%length))
|
2020-08-15 19:32:10 +05:30
|
|
|
allocate(source_thermal_externalheat_offset (phases%length), source=0)
|
|
|
|
|
2021-01-27 04:26:20 +05:30
|
|
|
do ph = 1, phases%length
|
|
|
|
phase => phases%get(ph)
|
|
|
|
if(count(mySources(:,ph)) == 0) cycle
|
2021-01-08 02:45:18 +05:30
|
|
|
thermal => phase%get('thermal')
|
|
|
|
sources => thermal%get('source')
|
2021-01-27 04:14:11 +05:30
|
|
|
do so = 1, sources%length
|
2021-01-27 04:26:20 +05:30
|
|
|
if(mySources(so,ph)) then
|
|
|
|
source_thermal_externalheat_offset(ph) = so
|
|
|
|
associate(prm => param(ph))
|
|
|
|
src => sources%get(so)
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-27 04:26:20 +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
|
|
|
|
2021-01-27 04:26:20 +05:30
|
|
|
prm%f_T = src%get_asFloats('f_T',requiredSize = size(prm%t_n))
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-27 04:26:20 +05:30
|
|
|
Nconstituents = count(material_phaseAt2 == ph)
|
2021-02-09 03:51:53 +05:30
|
|
|
call phase_allocateState(thermalState(ph)%p(so),Nconstituents,1,1,0)
|
2020-08-15 19:32:10 +05:30
|
|
|
end associate
|
|
|
|
endif
|
|
|
|
enddo
|
2019-02-23 01:07:41 +05:30
|
|
|
enddo
|
|
|
|
|
2021-01-26 12:25:06 +05:30
|
|
|
end function externalheat_init
|
2015-07-27 16:39:37 +05:30
|
|
|
|
2019-02-22 19:51:48 +05:30
|
|
|
|
2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2016-10-19 01:53:52 +05:30
|
|
|
!> @brief rate of change of state
|
|
|
|
!> @details state only contains current time to linearly interpolate given heat powers
|
2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-01-26 12:25:06 +05:30
|
|
|
module subroutine externalheat_dotState(ph, me)
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2019-03-24 15:48:59 +05:30
|
|
|
integer, intent(in) :: &
|
2021-01-19 15:02:56 +05:30
|
|
|
ph, &
|
|
|
|
me
|
2020-02-29 12:28:33 +05:30
|
|
|
|
2019-03-24 15:48:59 +05:30
|
|
|
integer :: &
|
2021-01-27 04:14:11 +05:30
|
|
|
so
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-27 04:14:11 +05:30
|
|
|
so = source_thermal_externalheat_offset(ph)
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-27 04:14:11 +05:30
|
|
|
thermalState(ph)%p(so)%dotState(1,me) = 1.0_pReal ! state is current time
|
2015-07-27 16:39:37 +05:30
|
|
|
|
2021-01-26 12:25:06 +05:30
|
|
|
end subroutine externalheat_dotState
|
2015-07-27 16:39:37 +05:30
|
|
|
|
2019-12-05 15:50:05 +05:30
|
|
|
|
2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-02-26 23:07:17 +05:30
|
|
|
!> @brief returns local heat generation rate
|
2015-07-27 16:39:37 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2021-01-27 04:26:20 +05:30
|
|
|
module subroutine externalheat_getRate(TDot, ph, me)
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2019-03-24 15:48:59 +05:30
|
|
|
integer, intent(in) :: &
|
2021-01-19 15:02:56 +05:30
|
|
|
ph, &
|
|
|
|
me
|
2019-03-24 15:48:59 +05:30
|
|
|
real(pReal), intent(out) :: &
|
2021-01-08 13:27:30 +05:30
|
|
|
TDot
|
2020-02-29 12:28:33 +05:30
|
|
|
|
2019-03-24 15:48:59 +05:30
|
|
|
integer :: &
|
2021-01-27 04:14:11 +05:30
|
|
|
so, interval
|
2019-03-24 15:48:59 +05:30
|
|
|
real(pReal) :: &
|
2020-02-26 23:07:17 +05:30
|
|
|
frac_time
|
|
|
|
|
2021-01-27 04:14:11 +05:30
|
|
|
so = source_thermal_externalheat_offset(ph)
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-27 04:26:20 +05:30
|
|
|
associate(prm => param(ph))
|
|
|
|
do interval = 1, prm%nIntervals ! scan through all rate segments
|
|
|
|
frac_time = (thermalState(ph)%p(so)%state(1,me) - prm%t_n(interval)) &
|
|
|
|
/ (prm%t_n(interval+1) - prm%t_n(interval)) ! fractional time within segment
|
|
|
|
if ( (frac_time < 0.0_pReal .and. interval == 1) &
|
|
|
|
.or. (frac_time >= 1.0_pReal .and. interval == prm%nIntervals) &
|
|
|
|
.or. (frac_time >= 0.0_pReal .and. frac_time < 1.0_pReal) ) &
|
|
|
|
TDot = prm%f_T(interval ) * (1.0_pReal - frac_time) + &
|
|
|
|
prm%f_T(interval+1) * frac_time ! interpolate heat rate between segment boundaries...
|
|
|
|
! ...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-27 04:26:20 +05:30
|
|
|
end subroutine externalheat_getRate
|
2020-02-26 23:07:17 +05:30
|
|
|
|
2021-01-26 05:50:45 +05:30
|
|
|
end submodule externalheat
|