wideBandAbsorptionEmission.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-2018 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 "basicSpecieMixture.H"
32#include "unitConversion.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 namespace radiation
39 {
41
43 (
47 );
48 }
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict,
57 const fvMesh& mesh
58)
59:
61 coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
62 speciesNames_(0),
63 specieIndex_(Zero),
64 lookUpTablePtr_(),
65 thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
66 Yj_(nSpecies_),
67 totalWaveLength_(0)
68{
69 label nBand = 0;
70 const dictionary& functionDicts = dict.optionalSubDict(typeName +"Coeffs");
71 for (const entry& dEntry : functionDicts)
72 {
73 if (!dEntry.isDict()) // safety
74 {
75 continue;
76 }
77
78 const dictionary& dict = dEntry.dict();
79
80 dict.readEntry("bandLimits", iBands_[nBand]);
81 dict.readEntry("EhrrCoeff", iEhrrCoeffs_[nBand]);
82 totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
83
84 label nSpec = 0;
85
86 const dictionary& specDicts = dict.subDict("species");
87 for (const entry& dEntry : specDicts)
88 {
89 const word& key = dEntry.keyword();
90
91 if (nBand == 0)
92 {
93 speciesNames_.insert(key, nSpec);
94 }
95 else if (!speciesNames_.found(key))
96 {
98 << "specie: " << key << " is not in all the bands"
99 << nl << exit(FatalError);
100 }
101 coeffs_[nBand][nSpec].initialise(specDicts.subDict(key));
102 nSpec++;
103 }
104 nBand++;
105 }
106 nBands_ = nBand;
107
108 if
109 (
110 coeffsDict_.found("lookUpTableFileName")
111 && "none" != coeffsDict_.get<word>("lookUpTableFileName")
112 )
113 {
114 lookUpTablePtr_.reset
115 (
117 (
118 coeffsDict_.get<fileName>("lookUpTableFileName"),
119 mesh.time().constant(),
120 mesh
121 )
122 );
123
124 if (!mesh.foundObject<volScalarField>("ft"))
125 {
127 << "specie ft is not present to use with "
128 << "lookUpTableFileName " << nl
129 << exit(FatalError);
130 }
131 }
132
133 // Check that all the species on the dictionary are present in the
134 // look-up table and save the corresponding indices of the look-up table
135
136 label j = 0;
137 forAllConstIters(speciesNames_, iter)
138 {
139 const word& specieName = iter.key();
140 const label index = iter.val();
141
142 volScalarField* fldPtr = mesh.getObjectPtr<volScalarField>(specieName);
143
144 if (lookUpTablePtr_)
145 {
146 if (lookUpTablePtr_().found(specieName))
147 {
148 const label fieldIndex =
149 lookUpTablePtr_().findFieldIndex(specieName);
150
151 Info<< "specie: " << specieName << " found on look-up table "
152 << " with index: " << fieldIndex << endl;
153
154 specieIndex_[index] = fieldIndex;
155 }
156 else if (fldPtr)
157 {
158 Yj_.set(j, fldPtr);
159 specieIndex_[index] = 0;
160 j++;
161 Info<< "specie: " << specieName << " is being solved" << endl;
162 }
163 else
164 {
166 << "specie: " << specieName
167 << " is neither in look-up table: "
168 << lookUpTablePtr_().tableName()
169 << " nor is being solved" << nl
170 << exit(FatalError);
171 }
172 }
173 else if (fldPtr)
174 {
175 Yj_.set(j, fldPtr);
176 specieIndex_[index] = 0;
177 j++;
178 }
179 else
180 {
182 << "There is no lookup table and the specie" << nl
183 << specieName << nl
184 << " is not found " << nl
185 << exit(FatalError);
186
187 }
188 }
189}
190
191
192// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
193
195{}
196
197
198// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199
202{
204 dynamic_cast<const basicSpecieMixture&>(thermo_);
205
206 const volScalarField& T = thermo_.T();
207 const volScalarField& p = thermo_.p();
208
210 (
212 (
214 (
215 "a",
216 mesh().time().timeName(),
217 mesh(),
220 ),
221 mesh(),
223 )
224 );
225
226 scalarField& a = ta.ref().primitiveFieldRef();
227
228 forAll(a, celli)
229 {
230 forAllConstIters(speciesNames_, iter)
231 {
232 const label n = iter();
233 scalar Xipi = 0;
234 if (specieIndex_[n] != 0)
235 {
236 const volScalarField& ft =
237 mesh_.lookupObject<volScalarField>("ft");
238
239 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
240
241 // moles*pressure [atm]
242 Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
243 }
244 else
245 {
246 scalar invWt = 0;
247 forAll(mixture.Y(), s)
248 {
249 invWt += mixture.Y(s)[celli]/mixture.W(s);
250 }
251
252 const label index = mixture.species().find(iter.key());
253
254 const scalar Xk =
255 mixture.Y(index)[celli]/(mixture.W(index)*invWt);
256
257 Xipi = Xk*paToAtm(p[celli]);
258 }
259
260 scalar Ti = T[celli];
261
263 coeffs_[bandi][n].coeffs(T[celli]);
264
265 if (coeffs_[bandi][n].invTemp())
266 {
267 Ti = 1.0/T[celli];
268 }
269
270 a[celli]+=
271 Xipi
272 *(
273 ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
274 + b[0]
275 );
276 }
277 }
278
279 return ta;
280}
281
282
285{
286 return aCont(bandi);
287}
288
289
292{
294 (
296 (
298 (
299 "E",
300 mesh().time().timeName(),
301 mesh(),
304 ),
305 mesh(),
307 )
308 );
309
310 const volScalarField* QdotPtr = mesh().findObject<volScalarField>("Qdot");
311
312 if (QdotPtr)
313 {
314 const volScalarField& Qdot = *QdotPtr;
315
316 if (Qdot.dimensions() == dimEnergy/dimTime)
317 {
318 E.ref().primitiveFieldRef() =
319 iEhrrCoeffs_[bandi]
320 *Qdot.primitiveField()
321 *(iBands_[bandi][1] - iBands_[bandi][0])
322 /totalWaveLength_
323 /mesh_.V();
324 }
325 else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
326 {
327 E.ref().primitiveFieldRef() =
328 iEhrrCoeffs_[bandi]
329 *Qdot.primitiveField()
330 *(iBands_[bandi][1] - iBands_[bandi][0])
331 /totalWaveLength_;
332 }
333 else
334 {
336 << "Incompatible dimensions for Qdot field" << endl;
337 }
338 }
339
340 return E;
341}
342
343
345(
348) const
349{
351
352 for (label j=0; j<nBands_; j++)
353 {
354 aLambda[j].primitiveFieldRef() = this->a(j);
355
356 a.primitiveFieldRef() +=
357 aLambda[j].primitiveField()
358 *(iBands_[j][1] - iBands_[j][0])
359 /totalWaveLength_;
360 }
361
362}
363
364
365// ************************************************************************* //
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.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
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:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
A class for handling file names.
Definition: fileName.H:76
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A list of lists. Interpolates based on the first dimension. The values must be positive and monotonic...
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Type * getObjectPtr(const word &name, const bool recursive=false) const
Model to supply absorption and emission coefficients for radiation modelling.
const fvMesh & mesh() const
Reference to the mesh.
const dictionary & dict() const
Reference to the dictionary.
wideBandAbsorptionEmission radiation absorption and emission coefficients for continuous phase.
FixedList< FixedList< absorptionCoeffs, nSpecies_ >, maxBands_ > coeffs_
Absorption coefficients.
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.
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const word dictName("faMeshDefinition")
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))
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
scalar Qdot
Definition: solveChemistry.H:2
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Unit conversion functions.