Source code for Mordicus.Modules.Safran.Containers.Loadings.Radiation
# -*- coding: utf-8 -*-
import os
from mpi4py import MPI
if MPI.COMM_WORLD.Get_size() > 1: # pragma: no cover
os.environ["OMP_NUM_THREADS"] = "1"
os.environ["OPENBLAS_NUM_THREADS"] = "1"
os.environ["MKL_NUM_THREADS"] = "1"
os.environ["VECLIB_MAXIMUM_THREADS"] = "1"
os.environ["NUMEXPR_NUM_THREADS"] = "1"
import numpy as np
from Mordicus.Core.Containers.Loadings.LoadingBase import LoadingBase
from Mordicus.Core.BasicAlgorithms import Interpolation as TI
[docs]class Radiation(LoadingBase):
"""
Class containing a Loading of type radiation boundary condition for thermal
problems.
Attributes
----------
TextTimes : np.ndarray or list of floats
time values on which external temperature values are provided
TextValues : np.ndarray or list of floats
external temperature values at the corresponding time values
StefanBoltzmannConstant : float
Stefan-Boltzmann constant used for the boundary condition computation
reducedPhiT : numpy.ndarray
of size (numberOfModes,) containing the precomputed vector
"red" for deriving the convection heat flux reduced external forces
contribution is the form: stefanBoltzmannConstant*(Text**4)*red
"""
def __init__(self, solutionName, set):
assert isinstance(set, str)
assert isinstance(solutionName, str)
assert solutionName == "T", "Radiation loading can only be applied on T solution types"
super(Radiation, self).__init__("T", set, "radiation")
self.TextTimes = None
self.TextValues = None
self.stefanBoltzmannConstant = 0.
self.reducedPhiT = None
[docs] def SetText(self, Text):
"""
Sets TextTimes and TextValues
Parameters
----------
Text : dict
dictionary with time steps (float) as keys and the values of the
external temperature values (float)
"""
# assert type of Text
assert isinstance(Text, dict)
assert np.all(
[isinstance(key, (float, np.float64)) for key in list(Text.keys())])
assert np.all([isinstance(key, (float, np.float64))
for key in list(Text.values())])
self.TextTimes = np.array(list(Text.keys()), dtype = float)
self.TextValues = np.array(list(Text.values()), dtype = float)
[docs] def SetStefanBoltzmannConstant(self, stefanBoltzmannConstant):
"""
Sets stefanBoltzmannConstant
Parameters
----------
stefanBoltzmannConstant : float
Stefan-Boltzmann constant used for the boundary condition computation
"""
# assert type of stefanBoltzmannConstant
assert isinstance(stefanBoltzmannConstant, (float, np.float64))
self.stefanBoltzmannConstant = stefanBoltzmannConstant
[docs] def GetTextAtTime(self, time: float)-> float:
"""
Computes and return Text at time, using PieceWiseLinearInterpolation
Parameters
----------
time : float
Returns
-------
float
Text at time
"""
# compute coefficient at time
Text = TI.PieceWiseLinearInterpolation(time,
self.TextTimes, self.TextValues)
return Text
[docs] def ReduceLoading(self, mesh, problemData, reducedOrderBases, operatorCompressionData):
"""
Computes and sets the reduced representation of the loading
Parameters
----------
mesh : BasicTools.Containers.UnstructuredMesh
mesh of the high-fidelity model
problemData : ProblemData
problemData containing the loading
reducedOrderBases : dict(str: np.ndarray)
dictionary with solutionNames (str) as keys and reducedOrderBases
(np.ndarray of size (numberOfModes, numberOfDOFs)) as values
operatorCompressionData : dict(str: custom_data_structure)
not used in this loading
dictionary with solutionNames (str) as keys and data structure
generated by the operator compression step as values
"""
from Mordicus.Modules.Safran.FE import FETools as FT
integrationWeights, phiAtIntegPoint = FT.ComputePhiAtIntegPoint(mesh, [self.GetSet()], relativeDimension = -1)
reducedPhiTAtIntegPoints = phiAtIntegPoint.dot(reducedOrderBases[self.solutionName].T)
self.reducedPhiT = np.einsum('tk,t->k', reducedPhiTAtIntegPoints, integrationWeights, optimize = True)
[docs] def ComputeContributionToReducedExternalForces(self, time):
"""
Computes and returns the reduced external forces contribution of the loading
Parameters
----------
time : float
Returns
-------
np.ndarray
of size (numberOfModes,)
"""
# assert type of time
assert isinstance(time, (float, np.float64))
Text = self.GetTextAtTime(time)
return self.stefanBoltzmannConstant*(Text**4)*self.reducedPhiT
def __str__(self):
res = "Radiation Loading with set "+self.GetSet()
return res
if __name__ == "__main__":# pragma: no cover
from Mordicus import RunTestFile
RunTestFile(__file__)