Loading...
Searching...
No Matches
greyMeanAbsorptionEmission.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) 2020 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"
33#include "basicSpecieMixture.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37namespace Foam
38{
39 namespace radiation
40 {
42
44 (
48 );
49 }
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const dictionary& dict,
58 const fvMesh& mesh
59)
60:
62 coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
63 speciesNames_(0),
64 specieIndex_(Zero),
65 lookUpTablePtr_(),
66 thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
67 EhrrCoeff_(coeffsDict_.get<scalar>("EhrrCoeff")),
68 Yj_(nSpecies_)
69{
70 if (!isA<basicSpecieMixture>(thermo_))
71 {
73 << "Model requires a multi-component thermo package"
74 << abort(FatalError);
75 }
76
77
78 label nFunc = 0;
79 const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
80
81 for (const entry& dEntry : functionDicts)
82 {
83 if (!dEntry.isDict()) // safety
84 {
85 continue;
86 }
87
88 const word& key = dEntry.keyword();
89 const dictionary& dict = dEntry.dict();
90
91 speciesNames_.insert(key, nFunc);
92
93 coeffs_[nFunc].initialise(dict);
94 nFunc++;
95 }
96
97 if
98 (
99 coeffsDict_.found("lookUpTableFileName")
100 && "none" != coeffsDict_.get<word>("lookUpTableFileName")
101 )
102 {
103 lookUpTablePtr_.reset
104 (
105 new interpolationLookUpTable<scalar>
106 (
107 coeffsDict_.get<fileName>("lookUpTableFileName"),
108 mesh.time().constant(),
109 mesh
110 )
111 );
112
113 if (!mesh.foundObject<volScalarField>("ft"))
114 {
115 FatalErrorInFunction
116 << "specie ft is not present to use with "
117 << "lookUpTableFileName " << nl
118 << exit(FatalError);
119 }
120 }
121
122 // Check that all the species on the dictionary are present in the
123 // look-up table and save the corresponding indices of the look-up table
124
125 label j = 0;
126 forAllConstIters(speciesNames_, iter)
127 {
128 const word& specieName = iter.key();
129 const label index = iter.val();
130
131 volScalarField* fldPtr = mesh.getObjectPtr<volScalarField>(specieName);
132
133 if (lookUpTablePtr_)
134 {
135 if (lookUpTablePtr_().found(specieName))
136 {
137 const label fieldIndex =
138 lookUpTablePtr_().findFieldIndex(specieName);
139
140 Info<< "specie: " << specieName << " found on look-up table "
141 << " with index: " << fieldIndex << endl;
142
143 specieIndex_[index] = fieldIndex;
144 }
145 else if (fldPtr)
146 {
147 Yj_.set(j, fldPtr);
148 specieIndex_[index] = 0;
149 j++;
150 Info<< "specie: " << iter.key() << " is being solved" << endl;
151 }
152 else
153 {
155 << "specie: " << specieName
156 << " is neither in look-up table: "
157 << lookUpTablePtr_().tableName()
158 << " nor is being solved" << nl
159 << exit(FatalError);
160 }
161 }
162 else if (fldPtr)
163 {
164 Yj_.set(j, fldPtr);
165 specieIndex_[index] = 0;
166 j++;
167 }
168 else
169 {
171 << "There is no lookup table and the specie" << nl
172 << specieName << nl
173 << " is not found " << nl
175 }
176 }
177}
178
179// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
180
182{}
183
184
185// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186
189{
191 dynamic_cast<const basicSpecieMixture&>(thermo_);
192
193 const volScalarField& T = thermo_.T();
194 const volScalarField& p = thermo_.p();
195
196
197 auto ta = volScalarField::New
198 (
199 "aCont" + name(bandI),
201 mesh(),
204 );
205 auto& a = ta.ref().primitiveFieldRef();
206
207 forAll(a, celli)
208 {
209 forAllConstIters(speciesNames_, iter)
210 {
211 label n = iter();
212 scalar Xipi = 0.0;
213 if (specieIndex_[n] != 0)
214 {
215 //Specie found in the lookUpTable.
216 const volScalarField& ft =
217 mesh_.lookupObject<volScalarField>("ft");
218
219 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
220 //moles x pressure [atm]
221 Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
222 }
223 else
224 {
225 scalar invWt = 0.0;
226 forAll(mixture.Y(), s)
227 {
228 invWt += mixture.Y(s)[celli]/mixture.W(s);
229 }
230
231 const label index = mixture.species().find(iter.key());
232 scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
233
234 Xipi = Xk*paToAtm(p[celli]);
235 }
236
237 const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
238
239 scalar Ti = T[celli];
240 // negative temperature exponents
241 if (coeffs_[n].invTemp())
242 {
243 Ti = 1.0/T[celli];
244 }
245 a[celli] +=
246 Xipi
247 *(
248 ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
249 + b[0]
250 );
251 }
253 ta.ref().correctBoundaryConditions();
254 return ta;
255}
256
257
260{
261 return aCont(bandI);
262}
263
264
267{
268 auto E = volScalarField::New
269 (
270 "ECont" + name(bandI),
272 mesh_,
274 );
275
276 const volScalarField* QdotPtr = mesh_.findObject<volScalarField>("Qdot");
277
278 if (QdotPtr)
279 {
280 const volScalarField& Qdot = *QdotPtr;
281
282 if (Qdot.dimensions() == dimEnergy/dimTime)
283 {
284 E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot/mesh_.V();
285 }
286 else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
287 {
288 E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot;
289 }
290 else
291 {
292 if (debug)
293 {
295 << "Incompatible dimensions for Qdot field" << endl;
296 }
297 }
298 }
299 else
300 {
302 << "Qdot field not found in mesh" << endl;
303 }
304
305 return E;
306}
307
308
309// ************************************************************************* //
bool found
label n
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).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
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
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
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.
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.
greyMeanAbsorptionEmission radiation absorption and emission coefficients for continuous phase
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.
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
A class for managing temporary objects.
Definition tmp.H:75
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
const volScalarField & T
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
auto & name
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.
Namespace for handling debugging switches.
Definition debug.C:45
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 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).
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...
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)
Definition errorManip.H:125
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).
Definition Ostream.H:50
dictionary dict
scalar Qdot
volScalarField & b
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
#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.