avoid mappings in bottom end functions
This commit is contained in:
parent
e0fa3e0b26
commit
ced7da4d62
|
@ -1091,6 +1091,8 @@ function homogenization_postResults(ip,el)
|
|||
use mesh, only: &
|
||||
mesh_element
|
||||
use material, only: &
|
||||
thermalMapping, &
|
||||
thermal_typeInstance, &
|
||||
mappingHomogenization, &
|
||||
homogState, &
|
||||
thermalState, &
|
||||
|
@ -1153,7 +1155,7 @@ function homogenization_postResults(ip,el)
|
|||
+ hydrogenfluxState(mappingHomogenization(2,ip,el))%sizePostResults) :: &
|
||||
homogenization_postResults
|
||||
integer(pInt) :: &
|
||||
startPos, endPos
|
||||
startPos, endPos, homog
|
||||
|
||||
homogenization_postResults = 0.0_pReal
|
||||
|
||||
|
@ -1184,11 +1186,13 @@ function homogenization_postResults(ip,el)
|
|||
case (THERMAL_isothermal_ID) chosenThermal
|
||||
|
||||
case (THERMAL_adiabatic_ID) chosenThermal
|
||||
homog = mappingHomogenization(2,ip,el)
|
||||
homogenization_postResults(startPos:endPos) = &
|
||||
thermal_adiabatic_postResults(ip, el)
|
||||
thermal_adiabatic_postResults(homog,thermal_typeInstance(homog),thermalMapping(homog)%p(ip,el))
|
||||
case (THERMAL_conduction_ID) chosenThermal
|
||||
homog = mappingHomogenization(2,ip,el)
|
||||
homogenization_postResults(startPos:endPos) = &
|
||||
thermal_conduction_postResults(ip, el)
|
||||
thermal_conduction_postResults(homog,thermal_typeInstance(homog),thermalMapping(homog)%p(ip,el))
|
||||
end select chosenThermal
|
||||
|
||||
startPos = endPos + 1_pInt
|
||||
|
|
|
@ -10,8 +10,6 @@ module thermal_adiabatic
|
|||
|
||||
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
|
||||
|
@ -98,7 +96,6 @@ subroutine thermal_adiabatic_init(fileUnit)
|
|||
maxNinstance = int(count(thermal_type == THERMAL_adiabatic_ID),pInt)
|
||||
if (maxNinstance == 0_pInt) return
|
||||
|
||||
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))
|
||||
thermal_adiabatic_output = ''
|
||||
|
@ -157,14 +154,13 @@ subroutine thermal_adiabatic_init(fileUnit)
|
|||
|
||||
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
|
||||
|
||||
! allocate state arrays
|
||||
sizeState = 1_pInt
|
||||
thermalState(section)%sizeState = sizeState
|
||||
thermalState(section)%sizePostResults = thermal_adiabatic_sizePostResults(instance)
|
||||
thermalState(section)%sizePostResults = sum(thermal_adiabatic_sizePostResult(:,instance))
|
||||
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))
|
||||
|
@ -344,6 +340,7 @@ function thermal_adiabatic_getSpecificHeat(ip,el)
|
|||
|
||||
end function thermal_adiabatic_getSpecificHeat
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns homogenized mass density
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -381,42 +378,38 @@ function thermal_adiabatic_getMassDensity(ip,el)
|
|||
thermal_adiabatic_getMassDensity/real(homogenization_Ngrains(mesh_element(3,el)),pReal)
|
||||
|
||||
end function thermal_adiabatic_getMassDensity
|
||||
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief return array of thermal results
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
function thermal_adiabatic_postResults(ip,el)
|
||||
function thermal_adiabatic_postResults(homog,instance,of) result(postResults)
|
||||
use material, only: &
|
||||
mappingHomogenization, &
|
||||
thermal_typeInstance, &
|
||||
thermalMapping, &
|
||||
temperature
|
||||
|
||||
implicit none
|
||||
integer(pInt), intent(in) :: &
|
||||
ip, & !< integration point
|
||||
el !< element
|
||||
real(pReal), dimension(thermal_adiabatic_sizePostResults(thermal_typeInstance(mappingHomogenization(2,ip,el)))) :: &
|
||||
thermal_adiabatic_postResults
|
||||
homog, &
|
||||
instance, &
|
||||
of
|
||||
|
||||
real(pReal), dimension(sum(thermal_adiabatic_sizePostResult(:,instance))) :: &
|
||||
postResults
|
||||
|
||||
integer(pInt) :: &
|
||||
instance, homog, offset, o, c
|
||||
|
||||
homog = mappingHomogenization(2,ip,el)
|
||||
offset = thermalMapping(homog)%p(ip,el)
|
||||
instance = thermal_typeInstance(homog)
|
||||
o, c
|
||||
|
||||
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)
|
||||
thermal_adiabatic_postResults(c+1_pInt) = temperature(homog)%p(offset)
|
||||
postResults(c+1_pInt) = temperature(homog)%p(of)
|
||||
c = c + 1
|
||||
end select
|
||||
enddo
|
||||
|
||||
end function thermal_adiabatic_postResults
|
||||
|
||||
end module thermal_adiabatic
|
||||
|
|
|
@ -10,8 +10,6 @@ module thermal_conduction
|
|||
|
||||
implicit none
|
||||
private
|
||||
integer(pInt), dimension(:), allocatable, public, protected :: &
|
||||
thermal_conduction_sizePostResults !< cumulative size of post results
|
||||
|
||||
integer(pInt), dimension(:,:), allocatable, target, public :: &
|
||||
thermal_conduction_sizePostResult !< size of each post result output
|
||||
|
@ -99,7 +97,6 @@ subroutine thermal_conduction_init(fileUnit)
|
|||
maxNinstance = int(count(thermal_type == THERMAL_conduction_ID),pInt)
|
||||
if (maxNinstance == 0_pInt) return
|
||||
|
||||
allocate(thermal_conduction_sizePostResults(maxNinstance), source=0_pInt)
|
||||
allocate(thermal_conduction_sizePostResult (maxval(homogenization_Noutput),maxNinstance),source=0_pInt)
|
||||
allocate(thermal_conduction_output (maxval(homogenization_Noutput),maxNinstance))
|
||||
thermal_conduction_output = ''
|
||||
|
@ -144,42 +141,40 @@ subroutine thermal_conduction_init(fileUnit)
|
|||
enddo parsingFile
|
||||
|
||||
initializeInstances: do section = 1_pInt, size(thermal_type)
|
||||
if (thermal_type(section) == THERMAL_conduction_ID) then
|
||||
NofMyHomog=count(material_homog==section)
|
||||
instance = thermal_typeInstance(section)
|
||||
if (thermal_type(section) /= THERMAL_conduction_ID) cycle
|
||||
NofMyHomog=count(material_homog==section)
|
||||
instance = thermal_typeInstance(section)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! Determine size of postResults array
|
||||
outputsLoop: do o = 1_pInt,thermal_conduction_Noutput(instance)
|
||||
select case(thermal_conduction_outputID(o,instance))
|
||||
case(temperature_ID)
|
||||
mySize = 1_pInt
|
||||
end select
|
||||
outputsLoop: do o = 1_pInt,thermal_conduction_Noutput(instance)
|
||||
select case(thermal_conduction_outputID(o,instance))
|
||||
case(temperature_ID)
|
||||
mySize = 1_pInt
|
||||
end select
|
||||
|
||||
if (mySize > 0_pInt) then ! any meaningful output found
|
||||
thermal_conduction_sizePostResult(o,instance) = mySize
|
||||
thermal_conduction_sizePostResults(instance) = thermal_conduction_sizePostResults(instance) + mySize
|
||||
endif
|
||||
enddo outputsLoop
|
||||
if (mySize > 0_pInt) then ! any meaningful output found
|
||||
thermal_conduction_sizePostResult(o,instance) = mySize
|
||||
endif
|
||||
enddo outputsLoop
|
||||
|
||||
! allocate state arrays
|
||||
sizeState = 0_pInt
|
||||
thermalState(section)%sizeState = sizeState
|
||||
thermalState(section)%sizePostResults = thermal_conduction_sizePostResults(instance)
|
||||
allocate(thermalState(section)%state0 (sizeState,NofMyHomog))
|
||||
allocate(thermalState(section)%subState0(sizeState,NofMyHomog))
|
||||
allocate(thermalState(section)%state (sizeState,NofMyHomog))
|
||||
sizeState = 0_pInt
|
||||
thermalState(section)%sizeState = sizeState
|
||||
thermalState(section)%sizePostResults = sum(thermal_conduction_sizePostResult(:,instance))
|
||||
allocate(thermalState(section)%state0 (sizeState,NofMyHomog))
|
||||
allocate(thermalState(section)%subState0(sizeState,NofMyHomog))
|
||||
allocate(thermalState(section)%state (sizeState,NofMyHomog))
|
||||
|
||||
nullify(thermalMapping(section)%p)
|
||||
thermalMapping(section)%p => mappingHomogenization(1,:,:)
|
||||
deallocate(temperature (section)%p)
|
||||
allocate (temperature (section)%p(NofMyHomog), source=thermal_initialT(section))
|
||||
deallocate(temperatureRate(section)%p)
|
||||
allocate (temperatureRate(section)%p(NofMyHomog), source=0.0_pReal)
|
||||
nullify(thermalMapping(section)%p)
|
||||
thermalMapping(section)%p => mappingHomogenization(1,:,:)
|
||||
deallocate(temperature (section)%p)
|
||||
allocate (temperature (section)%p(NofMyHomog), source=thermal_initialT(section))
|
||||
deallocate(temperatureRate(section)%p)
|
||||
allocate (temperatureRate(section)%p(NofMyHomog), source=0.0_pReal)
|
||||
|
||||
endif
|
||||
|
||||
enddo initializeInstances
|
||||
|
||||
end subroutine thermal_conduction_init
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -261,6 +256,7 @@ subroutine thermal_conduction_getSourceAndItsTangent(Tdot, dTdot_dT, T, ip, el)
|
|||
|
||||
end subroutine thermal_conduction_getSourceAndItsTangent
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns homogenized thermal conductivity in reference configuration
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -298,7 +294,8 @@ function thermal_conduction_getConductivity33(ip,el)
|
|||
thermal_conduction_getConductivity33/real(homogenization_Ngrains(mesh_element(3,el)),pReal)
|
||||
|
||||
end function thermal_conduction_getConductivity33
|
||||
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief returns homogenized specific heat capacity
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -374,7 +371,8 @@ function thermal_conduction_getMassDensity(ip,el)
|
|||
thermal_conduction_getMassDensity/real(homogenization_Ngrains(mesh_element(3,el)),pReal)
|
||||
|
||||
end function thermal_conduction_getMassDensity
|
||||
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief updates thermal state with solution from heat conduction PDE
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -403,41 +401,37 @@ subroutine thermal_conduction_putTemperatureAndItsRate(T,Tdot,ip,el)
|
|||
|
||||
end subroutine thermal_conduction_putTemperatureAndItsRate
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief return array of thermal results
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
function thermal_conduction_postResults(ip,el)
|
||||
function thermal_conduction_postResults(homog,instance,of) result(postResults)
|
||||
use material, only: &
|
||||
mappingHomogenization, &
|
||||
thermal_typeInstance, &
|
||||
temperature, &
|
||||
thermalMapping
|
||||
temperature
|
||||
|
||||
implicit none
|
||||
integer(pInt), intent(in) :: &
|
||||
ip, & !< integration point
|
||||
el !< element
|
||||
real(pReal), dimension(thermal_conduction_sizePostResults(thermal_typeInstance(mappingHomogenization(2,ip,el)))) :: &
|
||||
thermal_conduction_postResults
|
||||
homog, &
|
||||
instance, &
|
||||
of
|
||||
|
||||
real(pReal), dimension(sum(thermal_conduction_sizePostResult(:,instance))) :: &
|
||||
postResults
|
||||
|
||||
integer(pInt) :: &
|
||||
instance, homog, offset, o, c
|
||||
|
||||
homog = mappingHomogenization(2,ip,el)
|
||||
offset = thermalMapping(homog)%p(ip,el)
|
||||
instance = thermal_typeInstance(homog)
|
||||
o, c
|
||||
|
||||
c = 0_pInt
|
||||
thermal_conduction_postResults = 0.0_pReal
|
||||
|
||||
do o = 1_pInt,thermal_conduction_Noutput(instance)
|
||||
select case(thermal_conduction_outputID(o,instance))
|
||||
|
||||
case (temperature_ID)
|
||||
thermal_conduction_postResults(c+1_pInt) = temperature(homog)%p(offset)
|
||||
postResults(c+1_pInt) = temperature(homog)%p(of)
|
||||
c = c + 1
|
||||
end select
|
||||
enddo
|
||||
|
||||
end function thermal_conduction_postResults
|
||||
|
||||
end module thermal_conduction
|
||||
|
|
|
@ -26,14 +26,14 @@ subroutine thermal_isothermal_init()
|
|||
pInt
|
||||
use IO, only: &
|
||||
IO_timeStamp
|
||||
use config, only: &
|
||||
material_Nhomogenization
|
||||
use material
|
||||
use config
|
||||
|
||||
implicit none
|
||||
integer(pInt) :: &
|
||||
homog, &
|
||||
NofMyHomog, &
|
||||
sizeState
|
||||
NofMyHomog
|
||||
|
||||
write(6,'(/,a)') ' <<<+- thermal_'//THERMAL_isothermal_label//' init -+>>>'
|
||||
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
|
||||
|
@ -41,21 +41,19 @@ subroutine thermal_isothermal_init()
|
|||
|
||||
initializeInstances: do homog = 1_pInt, material_Nhomogenization
|
||||
|
||||
myhomog: if (thermal_type(homog) == THERMAL_isothermal_ID) then
|
||||
NofMyHomog = count(material_homog == homog)
|
||||
sizeState = 0_pInt
|
||||
thermalState(homog)%sizeState = sizeState
|
||||
thermalState(homog)%sizePostResults = sizeState
|
||||
allocate(thermalState(homog)%state0 (sizeState,NofMyHomog), source=0.0_pReal)
|
||||
allocate(thermalState(homog)%subState0(sizeState,NofMyHomog), source=0.0_pReal)
|
||||
allocate(thermalState(homog)%state (sizeState,NofMyHomog), source=0.0_pReal)
|
||||
if (thermal_type(homog) /= THERMAL_isothermal_ID) cycle
|
||||
NofMyHomog = count(material_homog == homog)
|
||||
thermalState(homog)%sizeState = 0_pInt
|
||||
thermalState(homog)%sizePostResults = 0_pInt
|
||||
allocate(thermalState(homog)%state0 (0_pInt,NofMyHomog), source=0.0_pReal)
|
||||
allocate(thermalState(homog)%subState0(0_pInt,NofMyHomog), source=0.0_pReal)
|
||||
allocate(thermalState(homog)%state (0_pInt,NofMyHomog), source=0.0_pReal)
|
||||
|
||||
deallocate(temperature (homog)%p)
|
||||
allocate (temperature (homog)%p(1), source=thermal_initialT(homog))
|
||||
deallocate(temperatureRate(homog)%p)
|
||||
allocate (temperatureRate(homog)%p(1), source=0.0_pReal)
|
||||
deallocate(temperature (homog)%p)
|
||||
allocate (temperature (homog)%p(1), source=thermal_initialT(homog))
|
||||
deallocate(temperatureRate(homog)%p)
|
||||
allocate (temperatureRate(homog)%p(1), source=0.0_pReal)
|
||||
|
||||
endif myhomog
|
||||
enddo initializeInstances
|
||||
|
||||
|
||||
|
|
Loading…
Reference in New Issue