Go to the documentation of this file.
45 #ifndef radiationModel_H
46 #define radiationModel_H
232 virtual bool read() = 0;
276 virtual label
nBands()
const = 0;
289 #define addToRadiationRunTimeSelectionTables(model) \
291 addToRunTimeSelectionTable \
298 addToRunTimeSelectionTable \
autoPtr< scatterModel > scatter_
Scatter model.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
label solverFreq_
Radiation solver frequency - number flow solver iterations per.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
dictionary coeffs_
Radiation model dictionary.
A class for managing temporary objects.
static const word relfectedFluxName_
Static name for reflected solar fluxes.
static const word primaryFluxName_
Static name for primary solar fluxes.
TypeName("radiationModel")
Runtime type information.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Abstract base-class for fluid and solid thermodynamic properties.
Forward declarations of fvMatrix specializations.
virtual void calculate()=0
Solve radiation equation(s)
const volScalarField & T_
Reference to the temperature field.
virtual bool read()=0
Read radiationProperties dictionary.
const sootModel & soot() const
Access to soot Model.
autoPtr< sootModel > soot_
Soot model.
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
declareRunTimeSelectionTable(autoPtr, radiationModel, T,(const volScalarField &T),(T))
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~radiationModel()
Destructor.
virtual tmp< volScalarField::Internal > Ru() const =0
Source term component (constant)
bool firstIter_
Flag to enable radiation model to be evaluated on first iteration.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
virtual label nBands() const =0
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Switch radiation() const
Radiation model on/off flag.
Mesh data needed to do the Finite Volume discretisation.
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
virtual void correct()
Main update/correction routine.
const fvMesh & mesh_
Reference to the mesh database.
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Macros to ease declaration of run-time selection tables.
Top level model for radiation modelling.
Base class for soor models.
virtual tmp< fvScalarMatrix > Sh(const basicThermo &thermo, const volScalarField &he) const
Energy source term.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
const Time & time_
Reference to the time database.
virtual tmp< volScalarField > Rp() const =0
Source term component (for power of T^4)
Model to supply absorption and emission coefficients for radiation modelling.
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
Switch radiation_
Radiation model on/off flag.
static const word externalRadHeatFieldName_
Static name external radiative fluxes.