2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2014-10-11 02:25:09 +05:30
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
2015-05-28 22:32:23 +05:30
!> @brief material subroutine for adiabatic temperature evolution
2014-10-09 19:36:45 +05:30
!> @details to be done
!--------------------------------------------------------------------------------------------------
module thermal_adiabatic
use prec , only : &
pReal , &
pInt
implicit none
private
integer ( pInt ) , dimension ( : ) , allocatable , public , protected :: &
thermal_adiabatic_sizePostResults !< cumulative size of post results
integer ( pInt ) , dimension ( : , : ) , allocatable , target , public :: &
thermal_adiabatic_sizePostResult !< size of each post result output
character ( len = 64 ) , dimension ( : , : ) , allocatable , target , public :: &
thermal_adiabatic_output !< name of each post result output
integer ( pInt ) , dimension ( : ) , allocatable , target , public :: &
2015-05-28 22:32:23 +05:30
thermal_adiabatic_Noutput !< number of outputs per instance of this thermal model
2014-10-09 19:36:45 +05:30
enum , bind ( c )
enumerator :: undefined_ID , &
temperature_ID
end enum
integer ( kind ( undefined_ID ) ) , dimension ( : , : ) , allocatable , private :: &
thermal_adiabatic_outputID !< ID of each post result output
public :: &
thermal_adiabatic_init , &
2015-05-28 22:32:23 +05:30
thermal_adiabatic_updateState , &
thermal_adiabatic_getSourceAndItsTangent , &
thermal_adiabatic_getSpecificHeat , &
thermal_adiabatic_getMassDensity , &
2014-10-09 19:36:45 +05:30
thermal_adiabatic_postResults
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2015-07-24 20:23:50 +05:30
subroutine thermal_adiabatic_init ( fileUnit )
2018-02-02 17:06:09 +05:30
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
2017-10-05 20:05:34 +05:30
use , intrinsic :: iso_fortran_env , only : &
compiler_version , &
compiler_options
#endif
2014-10-09 19:36:45 +05:30
use IO , only : &
IO_read , &
IO_lc , &
IO_getTag , &
IO_isBlank , &
IO_stringPos , &
IO_stringValue , &
IO_floatValue , &
IO_intValue , &
IO_warning , &
IO_error , &
IO_timeStamp , &
IO_EOF
2018-06-10 21:31:52 +05:30
use config_material , only : &
material_partHomogenization
2014-10-09 19:36:45 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
thermal_type , &
thermal_typeInstance , &
homogenization_Noutput , &
THERMAL_ADIABATIC_label , &
THERMAL_adiabatic_ID , &
material_homog , &
mappingHomogenization , &
2014-10-09 19:36:45 +05:30
thermalState , &
2015-05-28 22:32:23 +05:30
thermalMapping , &
2015-07-24 20:23:50 +05:30
thermal_initialT , &
2015-05-28 22:32:23 +05:30
temperature , &
2018-06-10 21:31:52 +05:30
temperatureRate
2014-10-09 19:36:45 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: fileUnit
2015-08-28 13:08:48 +05:30
integer ( pInt ) , allocatable , dimension ( : ) :: chunkPos
2015-05-28 22:32:23 +05:30
integer ( pInt ) :: maxNinstance , mySize = 0_pInt , section , instance , o
integer ( pInt ) :: sizeState
integer ( pInt ) :: NofMyHomog
2014-10-09 19:36:45 +05:30
character ( len = 65536 ) :: &
tag = '' , &
line = ''
2016-10-26 02:24:32 +05:30
write ( 6 , '(/,a)' ) ' <<<+- thermal_' / / THERMAL_ADIABATIC_label / / ' init -+>>>'
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2014-10-09 19:36:45 +05:30
#include "compilation_info.f90"
2015-05-28 22:32:23 +05:30
maxNinstance = int ( count ( thermal_type == THERMAL_adiabatic_ID ) , pInt )
2014-10-09 19:36:45 +05:30
if ( maxNinstance == 0_pInt ) return
2015-05-28 22:32:23 +05:30
allocate ( thermal_adiabatic_sizePostResults ( maxNinstance ) , source = 0_pInt )
allocate ( thermal_adiabatic_sizePostResult ( maxval ( homogenization_Noutput ) , maxNinstance ) , source = 0_pInt )
allocate ( thermal_adiabatic_output ( maxval ( homogenization_Noutput ) , maxNinstance ) )
2014-10-09 19:36:45 +05:30
thermal_adiabatic_output = ''
2015-05-28 22:32:23 +05:30
allocate ( thermal_adiabatic_outputID ( maxval ( homogenization_Noutput ) , maxNinstance ) , source = undefined_ID )
allocate ( thermal_adiabatic_Noutput ( maxNinstance ) , source = 0_pInt )
2014-10-09 19:36:45 +05:30
rewind ( fileUnit )
2015-05-28 22:32:23 +05:30
section = 0_pInt
do while ( trim ( line ) / = IO_EOF . and . IO_lc ( IO_getTag ( line , '<' , '>' ) ) / = material_partHomogenization ) ! wind forward to <homogenization>
2014-10-09 19:36:45 +05:30
line = IO_read ( fileUnit )
enddo
2015-05-28 22:32:23 +05:30
parsingFile : do while ( trim ( line ) / = IO_EOF ) ! read through sections of homog part
2014-10-09 19:36:45 +05:30
line = IO_read ( fileUnit )
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
if ( IO_getTag ( line , '<' , '>' ) / = '' ) then ! stop at next part
line = IO_read ( fileUnit , . true . ) ! reset IO_read
exit
endif
2015-05-28 22:32:23 +05:30
if ( IO_getTag ( line , '[' , ']' ) / = '' ) then ! next homog section
section = section + 1_pInt ! advance homog section counter
2014-10-09 19:36:45 +05:30
cycle ! skip to next line
endif
2015-05-28 22:32:23 +05:30
if ( section > 0_pInt ) then ; if ( thermal_type ( section ) == THERMAL_adiabatic_ID ) then ! do not short-circuit here (.and. with next if statemen). It's not safe in Fortran
2014-10-09 19:36:45 +05:30
2015-05-28 22:32:23 +05:30
instance = thermal_typeInstance ( section ) ! which instance of my thermal is present homog
2015-08-28 13:08:48 +05:30
chunkPos = IO_stringPos ( line )
tag = IO_lc ( IO_stringValue ( line , chunkPos , 1_pInt ) ) ! extract key
2014-10-09 19:36:45 +05:30
select case ( tag )
case ( '(output)' )
2015-08-28 13:08:48 +05:30
select case ( IO_lc ( IO_stringValue ( line , chunkPos , 2_pInt ) ) )
2014-10-09 19:36:45 +05:30
case ( 'temperature' )
thermal_adiabatic_Noutput ( instance ) = thermal_adiabatic_Noutput ( instance ) + 1_pInt
thermal_adiabatic_outputID ( thermal_adiabatic_Noutput ( instance ) , instance ) = temperature_ID
thermal_adiabatic_output ( thermal_adiabatic_Noutput ( instance ) , instance ) = &
2015-08-28 13:08:48 +05:30
IO_lc ( IO_stringValue ( line , chunkPos , 2_pInt ) )
2014-10-09 19:36:45 +05:30
end select
end select
endif ; endif
enddo parsingFile
2015-05-28 22:32:23 +05:30
initializeInstances : do section = 1_pInt , size ( thermal_type )
if ( thermal_type ( section ) == THERMAL_adiabatic_ID ) then
NofMyHomog = count ( material_homog == section )
instance = thermal_typeInstance ( section )
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
! Determine size of postResults array
outputsLoop : do o = 1_pInt , thermal_adiabatic_Noutput ( instance )
select case ( thermal_adiabatic_outputID ( o , instance ) )
case ( temperature_ID )
mySize = 1_pInt
end select
if ( mySize > 0_pInt ) then ! any meaningful output found
thermal_adiabatic_sizePostResult ( o , instance ) = mySize
thermal_adiabatic_sizePostResults ( instance ) = thermal_adiabatic_sizePostResults ( instance ) + mySize
endif
enddo outputsLoop
2015-05-28 22:32:23 +05:30
! allocate state arrays
sizeState = 1_pInt
thermalState ( section ) % sizeState = sizeState
thermalState ( section ) % sizePostResults = thermal_adiabatic_sizePostResults ( instance )
2015-07-24 20:23:50 +05:30
allocate ( thermalState ( section ) % state0 ( sizeState , NofMyHomog ) , source = thermal_initialT ( section ) )
allocate ( thermalState ( section ) % subState0 ( sizeState , NofMyHomog ) , source = thermal_initialT ( section ) )
allocate ( thermalState ( section ) % state ( sizeState , NofMyHomog ) , source = thermal_initialT ( section ) )
2015-05-28 22:32:23 +05:30
nullify ( thermalMapping ( section ) % p )
thermalMapping ( section ) % p = > mappingHomogenization ( 1 , : , : )
deallocate ( temperature ( section ) % p )
temperature ( section ) % p = > thermalState ( section ) % state ( 1 , : )
deallocate ( temperatureRate ( section ) % p )
allocate ( temperatureRate ( section ) % p ( NofMyHomog ) , source = 0.0_pReal )
2014-10-09 19:36:45 +05:30
endif
enddo initializeInstances
end subroutine thermal_adiabatic_init
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
!> @brief calculates adiabatic change in temperature based on local heat generation model
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
function thermal_adiabatic_updateState ( subdt , ip , el )
use numerics , only : &
err_thermal_tolAbs , &
err_thermal_tolRel
2015-03-18 23:33:18 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
mappingHomogenization , &
thermalState , &
temperature , &
temperatureRate , &
thermalMapping
2015-03-18 23:33:18 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
ip , & !< integration point number
el !< element number
real ( pReal ) , intent ( in ) :: &
subdt
2015-05-28 22:32:23 +05:30
logical , dimension ( 2 ) :: &
thermal_adiabatic_updateState
2015-03-18 23:33:18 +05:30
integer ( pInt ) :: &
2015-05-28 22:32:23 +05:30
homog , &
offset
real ( pReal ) :: &
T , Tdot , dTdot_dT
2015-03-18 23:33:18 +05:30
2015-05-28 22:32:23 +05:30
homog = mappingHomogenization ( 2 , ip , el )
offset = mappingHomogenization ( 1 , ip , el )
T = thermalState ( homog ) % subState0 ( 1 , offset )
call thermal_adiabatic_getSourceAndItsTangent ( Tdot , dTdot_dT , T , ip , el )
2015-07-27 16:39:37 +05:30
T = T + subdt * Tdot / ( thermal_adiabatic_getSpecificHeat ( ip , el ) * thermal_adiabatic_getMassDensity ( ip , el ) )
2015-03-18 23:33:18 +05:30
2015-05-28 22:32:23 +05:30
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 . ]
temperature ( homog ) % p ( thermalMapping ( homog ) % p ( ip , el ) ) = T
2015-08-05 16:39:38 +05:30
temperatureRate ( homog ) % p ( thermalMapping ( homog ) % p ( ip , el ) ) = &
( thermalState ( homog ) % state ( 1 , offset ) - thermalState ( homog ) % subState0 ( 1 , offset ) ) / ( subdt + tiny ( 0.0_pReal ) )
2015-05-28 22:32:23 +05:30
end function thermal_adiabatic_updateState
2015-03-18 23:33:18 +05:30
2014-11-01 00:33:08 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
!> @brief returns heat generation rate
2014-11-01 00:33:08 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
subroutine thermal_adiabatic_getSourceAndItsTangent ( Tdot , dTdot_dT , T , ip , el )
2014-11-01 00:33:08 +05:30
use math , only : &
2014-11-12 04:44:45 +05:30
math_Mandel6to33
2014-10-09 19:36:45 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
homogenization_Ngrains , &
2015-04-21 21:41:30 +05:30
mappingHomogenization , &
2016-05-27 15:16:34 +05:30
phaseAt , &
2015-05-28 22:32:23 +05:30
thermal_typeInstance , &
phase_Nsources , &
phase_source , &
2015-07-27 16:39:37 +05:30
SOURCE_thermal_dissipation_ID , &
SOURCE_thermal_externalheat_ID
2015-05-28 22:32:23 +05:30
use source_thermal_dissipation , only : &
source_thermal_dissipation_getRateAndItsTangent
2015-07-27 16:39:37 +05:30
use source_thermal_externalheat , only : &
source_thermal_externalheat_getRateAndItsTangent
2015-05-28 22:32:23 +05:30
use crystallite , only : &
crystallite_Tstar_v , &
crystallite_Lp
2014-10-09 19:36:45 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
ip , & !< integration point number
el !< element number
2015-05-28 22:32:23 +05:30
real ( pReal ) , intent ( in ) :: &
T
real ( pReal ) , intent ( out ) :: &
Tdot , dTdot_dT
real ( pReal ) :: &
my_Tdot , my_dTdot_dT
integer ( pInt ) :: &
phase , &
homog , &
offset , &
instance , &
grain , &
source
homog = mappingHomogenization ( 2 , ip , el )
offset = mappingHomogenization ( 1 , ip , el )
instance = thermal_typeInstance ( homog )
Tdot = 0.0_pReal
dTdot_dT = 0.0_pReal
do grain = 1 , homogenization_Ngrains ( homog )
2016-01-15 05:49:44 +05:30
phase = phaseAt ( grain , ip , el )
2015-05-28 22:32:23 +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_Tstar_v ( 1 : 6 , grain , ip , el ) , &
crystallite_Lp ( 1 : 3 , 1 : 3 , grain , ip , el ) , &
grain , ip , el )
2015-07-27 16:39:37 +05:30
case ( SOURCE_thermal_externalheat_ID )
call source_thermal_externalheat_getRateAndItsTangent ( my_Tdot , my_dTdot_dT , &
grain , ip , el )
2015-05-28 22:32:23 +05:30
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
2014-10-09 19:36:45 +05:30
2016-05-27 15:16:34 +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
2015-05-28 22:32:23 +05:30
end subroutine thermal_adiabatic_getSourceAndItsTangent
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
!> @brief returns homogenized specific heat capacity
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
function thermal_adiabatic_getSpecificHeat ( ip , el )
use lattice , only : &
lattice_specificHeat
2015-04-21 21:41:30 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
homogenization_Ngrains , &
mappingHomogenization , &
material_phase
use mesh , only : &
mesh_element
use crystallite , only : &
crystallite_push33ToRef
2015-04-21 21:41:30 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
ip , & !< integration point number
el !< element number
2015-05-28 22:32:23 +05:30
real ( pReal ) :: &
thermal_adiabatic_getSpecificHeat
integer ( pInt ) :: &
homog , grain
thermal_adiabatic_getSpecificHeat = 0.0_pReal
2015-04-21 21:41:30 +05:30
2015-05-28 22:32:23 +05:30
homog = mappingHomogenization ( 2 , ip , el )
do grain = 1 , homogenization_Ngrains ( mesh_element ( 3 , el ) )
thermal_adiabatic_getSpecificHeat = thermal_adiabatic_getSpecificHeat + &
lattice_specificHeat ( material_phase ( grain , ip , el ) )
enddo
thermal_adiabatic_getSpecificHeat = &
2016-05-27 15:16:34 +05:30
thermal_adiabatic_getSpecificHeat / real ( homogenization_Ngrains ( mesh_element ( 3 , el ) ) , pReal )
2015-04-21 21:41:30 +05:30
2015-05-28 22:32:23 +05:30
end function thermal_adiabatic_getSpecificHeat
2015-04-21 21:41:30 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
!> @brief returns homogenized mass density
2015-04-21 21:41:30 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
function thermal_adiabatic_getMassDensity ( ip , el )
use lattice , only : &
lattice_massDensity
2014-10-09 19:36:45 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
homogenization_Ngrains , &
mappingHomogenization , &
material_phase
use mesh , only : &
mesh_element
use crystallite , only : &
crystallite_push33ToRef
2014-10-09 19:36:45 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
ip , & !< integration point number
el !< element number
2015-05-28 22:32:23 +05:30
real ( pReal ) :: &
thermal_adiabatic_getMassDensity
integer ( pInt ) :: &
homog , grain
thermal_adiabatic_getMassDensity = 0.0_pReal
2014-10-09 19:36:45 +05:30
2015-05-28 22:32:23 +05:30
homog = mappingHomogenization ( 2 , ip , el )
2014-12-17 19:07:13 +05:30
2015-05-28 22:32:23 +05:30
do grain = 1 , homogenization_Ngrains ( mesh_element ( 3 , el ) )
thermal_adiabatic_getMassDensity = thermal_adiabatic_getMassDensity + &
lattice_massDensity ( material_phase ( grain , ip , el ) )
enddo
thermal_adiabatic_getMassDensity = &
2016-05-27 15:16:34 +05:30
thermal_adiabatic_getMassDensity / real ( homogenization_Ngrains ( mesh_element ( 3 , el ) ) , pReal )
2014-12-17 19:07:13 +05:30
2015-05-28 22:32:23 +05:30
end function thermal_adiabatic_getMassDensity
2014-12-17 19:07:13 +05:30
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
!> @brief return array of thermal results
2014-10-09 19:36:45 +05:30
!--------------------------------------------------------------------------------------------------
2015-05-28 22:32:23 +05:30
function thermal_adiabatic_postResults ( ip , el )
2014-10-09 19:36:45 +05:30
use material , only : &
2015-05-28 22:32:23 +05:30
mappingHomogenization , &
thermal_typeInstance , &
thermalMapping , &
temperature
2014-10-09 19:36:45 +05:30
implicit none
integer ( pInt ) , intent ( in ) :: &
ip , & !< integration point
el !< element
2015-05-28 22:32:23 +05:30
real ( pReal ) , dimension ( thermal_adiabatic_sizePostResults ( thermal_typeInstance ( mappingHomogenization ( 2 , ip , el ) ) ) ) :: &
2014-10-09 19:36:45 +05:30
thermal_adiabatic_postResults
integer ( pInt ) :: &
2015-05-28 22:32:23 +05:30
instance , homog , offset , o , c
2014-10-09 19:36:45 +05:30
2015-05-28 22:32:23 +05:30
homog = mappingHomogenization ( 2 , ip , el )
offset = thermalMapping ( homog ) % p ( ip , el )
instance = thermal_typeInstance ( homog )
2014-10-09 19:36:45 +05:30
c = 0_pInt
thermal_adiabatic_postResults = 0.0_pReal
do o = 1_pInt , thermal_adiabatic_Noutput ( instance )
select case ( thermal_adiabatic_outputID ( o , instance ) )
case ( temperature_ID )
2015-05-28 22:32:23 +05:30
thermal_adiabatic_postResults ( c + 1_pInt ) = temperature ( homog ) % p ( offset )
2014-10-09 19:36:45 +05:30
c = c + 1
end select
enddo
end function thermal_adiabatic_postResults
end module thermal_adiabatic