Go to the documentation of this file.
44 absorptionEmissionModel,
45 greyMeanAbsorptionEmission,
61 coeffsDict_((
dict.optionalSubDict(typeName +
"Coeffs"))),
66 EhrrCoeff_(coeffsDict_.get<scalar>(
"EhrrCoeff")),
69 if (!isA<basicSpecieMixture>(thermo_))
72 <<
"Model requires a multi-component thermo package"
78 const dictionary& functionDicts =
dict.optionalSubDict(typeName +
"Coeffs");
80 for (
const entry& dEntry : functionDicts)
87 const word& key = dEntry.keyword();
90 speciesNames_.insert(key, nFunc);
92 coeffs_[nFunc].initialise(
dict);
98 coeffsDict_.found(
"lookUpTableFileName")
99 &&
"none" != coeffsDict_.get<
word>(
"lookUpTableFileName")
102 lookUpTablePtr_.reset
106 coeffsDict_.get<
fileName>(
"lookUpTableFileName"),
107 mesh.time().constant(),
115 <<
"specie ft is not present to use with "
116 <<
"lookUpTableFileName " <<
nl
127 const word& specieName = iter.key();
128 const label index = iter.val();
132 if (!lookUpTablePtr_.empty())
134 if (lookUpTablePtr_().found(specieName))
136 const label fieldIndex =
137 lookUpTablePtr_().findFieldIndex(specieName);
139 Info<<
"specie: " << specieName <<
" found on look-up table "
140 <<
" with index: " << fieldIndex <<
endl;
142 specieIndex_[index] = fieldIndex;
147 specieIndex_[index] = 0;
149 Info<<
"specie: " << iter.key() <<
" is being solved" <<
endl;
154 <<
"specie: " << specieName
155 <<
" is neither in look-up table: "
156 << lookUpTablePtr_().tableName()
157 <<
" nor is being solved" <<
nl
164 specieIndex_[index] = 0;
170 <<
"There is no lookup table and the specie" <<
nl
172 <<
" is not found " <<
nl
202 "aCont" +
name(bandI),
210 extrapolatedCalculatedFvPatchVectorField::typeName
222 if (specieIndex_[
n] != 0)
228 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
230 Xipi = Ynft[specieIndex_[
n]]*
paToAtm(
p[celli]);
240 label index =
mixture.species()[iter.key()];
248 scalar Ti =
T[celli];
250 if (coeffs_[
n].invTemp())
257 ((((
b[5]*Ti +
b[4])*Ti +
b[3])*Ti +
b[2])*Ti +
b[1])*Ti
283 "ECont" +
name(bandI),
284 mesh_.time().timeName(),
302 E.
ref().primitiveFieldRef() = EhrrCoeff_*
Qdot/mesh_.V();
306 E.ref().primitiveFieldRef() = EhrrCoeff_*
Qdot;
313 <<
"Incompatible dimensions for Qdot field" <<
endl;
320 <<
"Qdot field not found in mesh" <<
endl;
int debug
Static debugging option.
A keyword and a list of tokens is an 'entry'.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
virtual ~greyMeanAbsorptionEmission()
Destructor.
A class for handling words, derived from Foam::string.
A class for handling file names.
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))
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimEnergy
Fundamental fluid thermodynamic properties.
const word dictName("blockMeshDict")
Unit conversion functions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
#define forAll(list, i)
Loop across all elements in list.
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
messageStream Info
Information stream (uses stdout - output is on the master only)
word name(const complex &c)
Return string representation of complex.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManip< error > abort(error &err)
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
forAllConstIters(mixture.phases(), phase)
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
Model to supply absorption and emission coefficients for radiation modelling.
const dimensionSet dimVolume(pow3(dimLength))
#define WarningInFunction
Report a warning using Foam::Warning.
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)