61 coeffsDict_((
dict.optionalSubDict(
typeName +
"Coeffs"))),
71 for (
const entry& dEntry : functionDicts)
80 dict.readEntry(
"bandLimits", iBands_[nBand]);
81 dict.readEntry(
"EhrrCoeff", iEhrrCoeffs_[nBand]);
82 totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
87 for (
const entry& dEntry : specDicts)
89 const word& key = dEntry.keyword();
93 speciesNames_.
insert(key, nSpec);
95 else if (!speciesNames_.found(key))
98 <<
"specie: " <<
key <<
" is not in all the bands"
101 coeffs_[nBand][nSpec].initialise(specDicts.subDict(key));
110 coeffsDict_.found(
"lookUpTableFileName")
111 &&
"none" != coeffsDict_.get<word>(
"lookUpTableFileName")
114 lookUpTablePtr_.reset
116 new interpolationLookUpTable<scalar>
118 coeffsDict_.get<fileName>(
"lookUpTableFileName"),
119 mesh.time().constant(),
124 if (!mesh.foundObject<volScalarField>(
"ft"))
127 <<
"specie ft is not present to use with "
128 <<
"lookUpTableFileName " << nl
139 const word& specieName = iter.key();
140 const label index = iter.val();
146 if (lookUpTablePtr_().
found(specieName))
148 const label fieldIndex =
149 lookUpTablePtr_().findFieldIndex(specieName);
151 Info<<
"specie: " << specieName <<
" found on look-up table "
152 <<
" with index: " << fieldIndex <<
endl;
154 specieIndex_[index] = fieldIndex;
159 specieIndex_[index] = 0;
161 Info<<
"specie: " << specieName <<
" is being solved" <<
endl;
166 <<
"specie: " << specieName
167 <<
" is neither in look-up table: "
168 << lookUpTablePtr_().tableName()
169 <<
" nor is being solved" <<
nl
176 specieIndex_[index] = 0;
182 <<
"There is no lookup table and the specie" <<
nl
184 <<
" is not found " <<
nl
222 const label
n = iter();
224 if (specieIndex_[
n] != 0)
229 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
232 Xipi = Ynft[specieIndex_[
n]]*
paToAtm(
p[celli]);
242 const label index =
mixture.species().find(iter.key());
250 scalar Ti =
T[celli];
253 coeffs_[bandi][
n].coeffs(
T[celli]);
255 if (coeffs_[bandi][
n].invTemp())
263 ((((
b[5]*Ti +
b[4])*Ti +
b[3])*Ti +
b[2])*Ti +
b[1])*Ti
299 E.ref().primitiveFieldRef() =
301 *
Qdot.primitiveField()
302 *(iBands_[bandi][1] - iBands_[bandi][0])
308 E.ref().primitiveFieldRef() =
310 *
Qdot.primitiveField()
311 *(iBands_[bandi][1] - iBands_[bandi][0])
317 <<
"Incompatible dimensions for Qdot field" <<
endl;
328 PtrList<volScalarField>& aLambda
333 for (label j=0; j<nBands_; j++)
335 aLambda[j].primitiveFieldRef() = this->a(j);
337 a.primitiveFieldRef() +=
338 aLambda[j].primitiveField()
339 *(iBands_[j][1] - iBands_[j][0])
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())
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
@ 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...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
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,...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
A keyword and a list of tokens is an 'entry'.
Fundamental fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
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.
wideBandAbsorptionEmission radiation absorption and emission coefficients for continuous phase.
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
FixedList< FixedList< absorptionCoeffs, nSpecies_ >, maxBands_ > coeffs_
Absorption coefficients.
tmp< volScalarField > aCont(const label bandi=0) const
Absorption coefficient for continuous phase.
wideBandAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
static const int nSpecies_
Maximum number of species considered for absorptivity.
tmp< volScalarField > ECont(const label bandi=0) const
Emission contribution for continuous phase.
virtual ~wideBandAbsorptionEmission()
Destructor.
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
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.
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).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
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.