35 namespace radiationModels
37 namespace absorptionEmissionModels
63 specieIndex_(
label(0)),
68 if (!isA<fluidMulticomponentThermo>(
thermo_))
71 <<
"Model requires a multi-component thermo package"
83 const word& key = iter().keyword();
90 if (
dict.found(
"lookUpTableFileName"))
108 <<
"specie ft is not present to use with "
109 <<
"lookUpTableFileName " <<
nl
127 Info<<
"specie: " << iter.key() <<
" found on look-up table "
128 <<
" with index: " << index <<
endl;
137 Info<<
"specie: " << iter.key() <<
" is being solved" <<
endl;
142 <<
"specie: " << iter.key()
143 <<
" is neither in look-up table: "
145 <<
" nor is being solved" <<
nl
158 <<
" there is not lookup table and the specie" <<
nl
160 <<
" is not found " <<
nl
192 "aCont" +
name(bandI),
195 extrapolatedCalculatedFvPatchVectorField::typeName
209 if (specieIndex_[
n] != 0)
215 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
218 Xipi = unitAtm.
toUser(Ynft[specieIndex_[
n]]*
p[celli]);
225 invWt += mcThermo.
Y(
s)[celli]/mcThermo.
WiValue(
s);
231 mcThermo.
Y(index)[celli]/(mcThermo.
WiValue(index)*invWt);
233 Xipi = unitAtm.
toUser(Xk*
p[celli]);
238 scalar Ti =
T[celli];
240 if (coeffs_[
n].invTemp())
247 ((((
b[5]*Ti +
b[4])*Ti +
b[3])*Ti +
b[2])*Ti +
b[1])*Ti
252 ta.
ref().correctBoundaryConditions();
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
static const word & constant()
Return constant name.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A class for handling file names.
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
virtual const speciesTable & species() const =0
The table of species.
virtual scalar WiValue(const label speciei) const =0
Molecular weight [kg/kmol].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
bool foundObject(const word &name) const
Is the named Type in registry.
A base class for physical properties.
void initialise(const dictionary &)
Model to supply absorption and emission coefficients for radiation modelling.
const fvMesh & mesh() const
Reference to the mesh.
virtual tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
const fluidThermo & thermo_
Thermo package.
UPtrList< volScalarField > Yj_
Pointer list of species in the registry involved in the absorption.
HashTable< label > speciesNames_
Hash table of species names.
autoPtr< interpolationLookUpTable > lookUpTablePtr_
Look-up table of species related to ft.
virtual ~greyMean()
Destructor.
FixedList< label, nSpecies_ > specieIndex_
Indices of species in the look-up table.
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.
greyMean(const dictionary &dict, const fvMesh &mesh, const word &modelName=typeName)
Construct from components.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
T toUser(const T &) const
Convert a value to user units.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
defineTypeNameAndDebug(greyMeanCombustion, 0)
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimLength
const HashTable< unitConversion > & units()
Get the table of unit conversions.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.