DAMASK_EICMD/src/phase_thermal_externalheat.f90

134 lines
4.9 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 :: &
2021-02-13 23:11:30 +05:30
source_thermal_externalheat_offset !< which source is my current thermal dissipation 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
2021-03-05 01:46:36 +05:30
integer :: so,Nmembers,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
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-03-11 22:30:07 +05:30
prm%t_n = src%get_as1dFloat('t_n')
2021-01-27 04:26:20 +05:30
prm%nIntervals = size(prm%t_n) - 1
2020-02-26 23:07:17 +05:30
2021-03-11 22:30:07 +05:30
prm%f_T = src%get_as1dFloat('f_T',requiredSize = size(prm%t_n))
2020-02-26 23:07:17 +05:30
2021-04-06 15:08:44 +05:30
Nmembers = count(material_phaseID == ph)
2021-03-05 01:46:36 +05:30
call phase_allocateState(thermalState(ph)%p(so),Nmembers,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-04-25 11:36:52 +05:30
module subroutine externalheat_dotState(ph, en)
2020-02-26 23:07:17 +05:30
integer, intent(in) :: &
2021-01-19 15:02:56 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
2020-02-29 12:28:33 +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-04-25 11:36:52 +05:30
thermalState(ph)%p(so)%dotState(1,en) = 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-04-25 11:36:52 +05:30
module function externalheat_f_T(ph,en) result(f_T)
2020-02-26 23:07:17 +05:30
integer, intent(in) :: &
2021-01-19 15:02:56 +05:30
ph, &
2021-04-25 11:36:52 +05:30
en
2021-04-08 02:11:49 +05:30
real(pReal) :: &
f_T
2020-02-29 12:28:33 +05:30
integer :: &
2021-01-27 04:14:11 +05:30
so, interval
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))
2021-03-25 23:52:59 +05:30
do interval = 1, prm%nIntervals ! scan through all rate segments
2021-04-25 11:36:52 +05:30
frac_time = (thermalState(ph)%p(so)%state(1,en) - prm%t_n(interval)) &
2021-03-25 23:52:59 +05:30
/ (prm%t_n(interval+1) - prm%t_n(interval)) ! fractional time within segment
2021-01-27 04:26:20 +05:30
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) ) &
2021-04-08 02:11:49 +05:30
f_T = prm%f_T(interval ) * (1.0_pReal - frac_time) + &
prm%f_T(interval+1) * frac_time ! interpolate heat rate between segment boundaries...
2021-04-25 11:36:52 +05:30
! ...or extrapolate if outside of bounds
2021-01-27 04:26:20 +05:30
enddo
2020-02-29 12:28:33 +05:30
end associate
2020-02-26 23:07:17 +05:30
2021-04-08 02:11:49 +05:30
end function externalheat_f_T
2020-02-26 23:07:17 +05:30
2021-01-26 05:50:45 +05:30
end submodule externalheat