Go to the documentation of this file.
44 absorptionEmissionModel,
45 greyMeanSolidAbsorptionEmission,
54 greyMeanSolidAbsorptionEmission::X(
const word specie)
const
72 Yi[iCell]/mixture_.
rho(specieI,
p[iCell],
T[iCell]);
76 const label mySpecieI = mixture_.
species()[specie];
79 Xj[iCell] = Yj[iCell]/mixture_.
rho(mySpecieI,
p[iCell],
T[iCell]);
95 coeffsDict_((
dict.optionalSubDict(typeName +
"Coeffs"))),
99 solidData_(mixture_.Y().size())
101 if (!isA<basicSpecieMixture>(thermo_))
104 <<
"Model requires a multi-component thermo package"
109 const dictionary& functionDicts =
dict.optionalSubDict(typeName +
"Coeffs");
111 for (
const entry& dEntry : functionDicts)
113 if (!dEntry.isDict())
118 const word&
key = dEntry.keyword();
121 if (!mixture_.contains(
key))
124 <<
" specie: " <<
key <<
" is not found in the solid mixture"
126 <<
" specie is the mixture are:" << mixture_.species() <<
nl
129 speciesNames_.insert(
key, nFunc);
131 dict.readEntry(
"absorptivity", solidData_[nFunc][absorptivity]);
132 dict.readEntry(
"emissivity", solidData_[nFunc][emissivity]);
142 Foam::radiation::greyMeanSolidAbsorptionEmission::
143 calc(
const label propertyId)
const
159 extrapolatedCalculatedFvPatchVectorField::typeName
167 if (mixture_.contains(iter.key()))
169 a += solidData_[iter()][propertyId]*X(iter.key());
173 ta.ref().correctBoundaryConditions();
184 return calc(emissivity);
194 return calc(absorptivity);
A keyword and a list of tokens is an 'entry'.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m3].
A class for handling words, derived from Foam::string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
virtual volScalarField & p()
Pressure [Pa].
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const speciesTable & species() const
Return the table of species.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Unit conversion functions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
#define forAll(list, i)
Loop across all elements in list.
Fundamental solid thermodynamic properties.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual const volScalarField & T() const
Temperature [K].
forAllConstIters(mixture.phases(), phase)
Model to supply absorption and emission coefficients for radiation modelling.
#define WarningInFunction
Report a warning using Foam::Warning.
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
const dimensionSet dimless
Dimensionless.
word dictName() const
The local dictionary name (final part of scoped name)