42 absorptionEmissionModel,
43 greyMeanSolidAbsorptionEmission,
52 greyMeanSolidAbsorptionEmission::X(
const word specie)
const 57 tmp<scalarField> tXj(
new scalarField(T.primitiveField().size(), 0.0));
60 tmp<scalarField> tRhoInv(
new scalarField(T.primitiveField().size(), 0.0));
63 forAll(mixture_.Y(), specieI)
70 Yi[iCell]/mixture_.rho(specieI, p[iCell], T[iCell]);
74 const label mySpecieI = mixture_.species()[specie];
77 Xj[iCell] = Yj[iCell]/mixture_.rho(mySpecieI, p[iCell], T[iCell]);
96 mixture_(dynamic_cast<const basicSpecieMixture&>(thermo_)),
97 solidData_(mixture_.Y().size())
99 if (!isA<basicSpecieMixture>(thermo_))
102 <<
"Model requires a multi-component thermo package" 112 if (!iter().isDict())
116 const word& key = iter().keyword();
117 if (!mixture_.contains(key))
120 <<
" specie: " << key <<
" is not found in the solid mixture" 122 <<
" specie is the mixture are:" << mixture_.species() <<
nl 125 speciesNames_.insert(key, nFunc);
127 dict.
lookup(
"absorptivity") >> solidData_[nFunc][absorptivity];
128 dict.
lookup(
"emissivity") >> solidData_[nFunc][emissivity];
144 Foam::radiation::greyMeanSolidAbsorptionEmission::
145 calc(
const label propertyId)
const 161 extrapolatedCalculatedFvPatchVectorField::typeName
169 if (mixture_.contains(iter.key()))
171 a += solidData_[iter()][propertyId]*X(iter.key());
175 ta.
ref().correctBoundaryConditions();
186 return calc(emissivity);
196 return calc(absorptivity);
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
T & ref() const
Return non-const reference or generate a fatal error.
Unit conversion functions.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Model to supply absorption and emission coefficients for radiation modelling.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
errorManip< error > abort(error &err)
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Fundamental solid thermodynamic properties.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
A class for managing temporary objects.
virtual ~greyMeanSolidAbsorptionEmission()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.