62 coeffsDict_((
dict.optionalSubDict(
typeName +
"Coeffs"))),
67 EhrrCoeff_(coeffsDict_.get<scalar>(
"EhrrCoeff")),
73 <<
"Model requires a multi-component thermo package"
81 for (
const entry& dEntry : functionDicts)
88 const word&
key = dEntry.keyword();
91 speciesNames_.insert(key, nFunc);
99 coeffsDict_.found(
"lookUpTableFileName")
100 &&
"none" != coeffsDict_.get<word>(
"lookUpTableFileName")
103 lookUpTablePtr_.reset
105 new interpolationLookUpTable<scalar>
107 coeffsDict_.get<fileName>(
"lookUpTableFileName"),
108 mesh.time().constant(),
113 if (!mesh.foundObject<volScalarField>(
"ft"))
116 <<
"specie ft is not present to use with "
117 <<
"lookUpTableFileName " << nl
128 const word& specieName = iter.key();
129 const label index = iter.val();
135 if (lookUpTablePtr_().
found(specieName))
137 const label fieldIndex =
138 lookUpTablePtr_().findFieldIndex(specieName);
140 Info<<
"specie: " << specieName <<
" found on look-up table "
141 <<
" with index: " << fieldIndex <<
endl;
143 specieIndex_[index] = fieldIndex;
148 specieIndex_[index] = 0;
150 Info<<
"specie: " << iter.key() <<
" is being solved" <<
endl;
155 <<
"specie: " << specieName
156 <<
" is neither in look-up table: "
157 << lookUpTablePtr_().tableName()
158 <<
" nor is being solved" <<
nl
165 specieIndex_[index] = 0;
171 <<
"There is no lookup table and the specie" <<
nl
173 <<
" is not found " <<
nl
199 "aCont" +
name(bandI),
205 auto& a = ta.ref().primitiveFieldRef();
213 if (specieIndex_[
n] != 0)
219 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
221 Xipi = Ynft[specieIndex_[
n]]*
paToAtm(
p[celli]);
231 const label index =
mixture.species().find(iter.key());
239 scalar Ti =
T[celli];
241 if (coeffs_[
n].invTemp())
248 ((((
b[5]*Ti +
b[4])*Ti +
b[3])*Ti +
b[2])*Ti +
b[1])*Ti
253 ta.ref().correctBoundaryConditions();
270 "ECont" +
name(bandI),
284 E.ref().primitiveFieldRef() = EhrrCoeff_*
Qdot/mesh_.V();
288 E.ref().primitiveFieldRef() = EhrrCoeff_*
Qdot;
295 <<
"Incompatible dimensions for Qdot field" <<
endl;
302 <<
"Qdot field not found in mesh" <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
label find(const T &val) const
Find index of the first occurrence of the value.
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Abstract base-class for fluid and solid thermodynamic properties.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A keyword and a list of tokens is an 'entry'.
Fundamental fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
FixedList< scalar, nCoeffs_ > coeffArray
Model to supply absorption and emission coefficients for radiation modelling.
const fvMesh & mesh_
Reference to the fvMesh.
virtual tmp< volScalarField > E(const label bandI=0) const
Emission contribution (net).
virtual tmp< volScalarField > a(const label bandI=0) const
Absorption coefficient (net).
const fvMesh & mesh() const
Reference to the mesh.
absorptionEmissionModel(const dictionary &dict, const fvMesh &mesh)
Construct from components.
const dictionary & dict() const
Reference to the dictionary.
greyMeanAbsorptionEmission radiation absorption and emission coefficients for continuous phase
static const int nSpecies_
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
absorptionCoeffs coeffs_[nSpecies_]
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
virtual ~greyMeanAbsorptionEmission()
Destructor.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const word dictName("faMeshDefinition")
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for radiation modelling.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Unit conversion functions.