Loading...
Searching...
No Matches
greyMeanSolidAbsorptionEmission.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "unitConversion.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36namespace Foam
37{
38 namespace radiation
39 {
41
43 (
47 );
48 }
49}
50
51// * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
52
53Foam::tmp<Foam::scalarField> Foam::radiation::
54greyMeanSolidAbsorptionEmission::X(const word specie) const
55{
56 const volScalarField& T = thermo_.T();
57 const volScalarField& p = thermo_.p();
58
59 auto tXj = tmp<scalarField>::New(T.primitiveField().size(), Zero);
60 auto& Xj = tXj.ref();
61
62 auto tRhoInv = tmp<scalarField>::New(T.primitiveField().size(), Zero);
63 auto& rhoInv = tRhoInv.ref();
64
65 forAll(mixture_.Y(), specieI)
66 {
67 const scalarField& Yi = mixture_.Y()[specieI];
68
69 forAll(rhoInv, iCell)
70 {
71 rhoInv[iCell] +=
72 Yi[iCell]/mixture_.rho(specieI, p[iCell], T[iCell]);
73 }
74 }
75 const scalarField& Yj = mixture_.Y(specie);
76 const label mySpecieI = mixture_.species().find(specie);
77 forAll(Xj, iCell)
78 {
79 Xj[iCell] = Yj[iCell]/mixture_.rho(mySpecieI, p[iCell], T[iCell]);
80 }
81
82 return (Xj/rhoInv);
83}
84
85// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86
89(
90 const dictionary& dict,
91 const fvMesh& mesh
92)
93:
95 coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
96 thermo_(mesh.lookupObject<solidThermo>(basicThermo::dictName)),
97 speciesNames_(0),
98 mixture_(dynamic_cast<const basicSpecieMixture&>(thermo_)),
99 solidData_(mixture_.Y().size())
100{
101 if (!isA<basicSpecieMixture>(thermo_))
102 {
104 << "Model requires a multi-component thermo package"
105 << abort(FatalError);
106 }
107
108 label nFunc = 0;
109 const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
110
111 for (const entry& dEntry : functionDicts)
112 {
113 if (!dEntry.isDict()) // safety
114 {
115 continue;
116 }
117
118 const word& key = dEntry.keyword();
119 const dictionary& dict = dEntry.dict();
120
121 if (!mixture_.contains(key))
122 {
124 << " specie: " << key << " is not found in the solid mixture"
125 << nl
126 << " specie is the mixture are:" << mixture_.species() << nl
127 << nl << endl;
128 }
129 speciesNames_.insert(key, nFunc);
130
131 dict.readEntry("absorptivity", solidData_[nFunc][absorptivity]);
132 dict.readEntry("emissivity", solidData_[nFunc][emissivity]);
133
134 nFunc++;
135 }
136}
137
138
139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140
141Foam::tmp<Foam::volScalarField>
142Foam::radiation::greyMeanSolidAbsorptionEmission::
143calc(const label propertyId) const
144{
145 auto ta = volScalarField::New
146 (
147 "a",
149 mesh(),
152 );
153 auto& a = ta.ref().primitiveFieldRef();
154
155 forAllConstIters(speciesNames_, iter)
156 {
157 if (mixture_.contains(iter.key()))
158 {
159 a += solidData_[iter.val()][propertyId]*X(iter.key());
160 }
161 }
163 ta.ref().correctBoundaryConditions();
164 return ta;
165}
166
167
170(
171 const label bandI
172) const
173{
174 return calc(emissivity);
175}
176
177
180(
181 const label bandI
182) const
183{
184 return calc(absorptivity);
185}
186
187
188// ************************************************************************* //
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).
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Model to supply absorption and emission coefficients for radiation modelling.
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.
greyMeanSolidAbsorptionEmission radiation absorption and emission coefficients for solid mixture
greyMeanSolidAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Fundamental solid thermodynamic properties.
Definition solidThermo.H:51
Base class of the thermophysical property types.
Definition specie.H:64
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
#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.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.