Loading...
Searching...
No Matches
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-2023 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"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
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 (
116 new interpolationLookUpTable<scalar>
117 (
118 coeffsDict_.get<fileName>("lookUpTableFileName"),
119 mesh.time().constant(),
120 mesh
121 )
122 );
123
124 if (!mesh.foundObject<volScalarField>("ft"))
125 {
126 FatalErrorInFunction
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
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
209 auto ta = volScalarField::New
210 (
211 "a",
213 mesh(),
215 );
216 scalarField& a = ta.ref().primitiveFieldRef();
217
218 forAll(a, celli)
219 {
220 forAllConstIters(speciesNames_, iter)
221 {
222 const label n = iter();
223 scalar Xipi = 0;
224 if (specieIndex_[n] != 0)
225 {
226 const volScalarField& ft =
227 mesh_.lookupObject<volScalarField>("ft");
228
229 const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
230
231 // moles*pressure [atm]
232 Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
233 }
234 else
235 {
236 scalar invWt = 0;
237 forAll(mixture.Y(), s)
238 {
239 invWt += mixture.Y(s)[celli]/mixture.W(s);
240 }
241
242 const label index = mixture.species().find(iter.key());
243
244 const scalar Xk =
245 mixture.Y(index)[celli]/(mixture.W(index)*invWt);
246
247 Xipi = Xk*paToAtm(p[celli]);
248 }
249
250 scalar Ti = T[celli];
251
253 coeffs_[bandi][n].coeffs(T[celli]);
254
255 if (coeffs_[bandi][n].invTemp())
256 {
257 Ti = 1.0/T[celli];
258 }
259
260 a[celli]+=
261 Xipi
262 *(
263 ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
264 + b[0]
265 );
266 }
268
269 return ta;
270}
271
272
275{
276 return aCont(bandi);
277}
278
279
282{
283 auto E = volScalarField::New
284 (
285 "E",
287 mesh(),
289 );
290
291 const volScalarField* QdotPtr = mesh().findObject<volScalarField>("Qdot");
292
293 if (QdotPtr)
294 {
295 const volScalarField& Qdot = *QdotPtr;
296
297 if (Qdot.dimensions() == dimEnergy/dimTime)
298 {
299 E.ref().primitiveFieldRef() =
300 iEhrrCoeffs_[bandi]
301 *Qdot.primitiveField()
302 *(iBands_[bandi][1] - iBands_[bandi][0])
303 /totalWaveLength_
304 /mesh_.V();
305 }
306 else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
307 {
308 E.ref().primitiveFieldRef() =
309 iEhrrCoeffs_[bandi]
310 *Qdot.primitiveField()
311 *(iBands_[bandi][1] - iBands_[bandi][0])
312 /totalWaveLength_;
313 }
314 else
315 {
317 << "Incompatible dimensions for Qdot field" << endl;
319 }
320
321 return E;
322}
323
324
326(
328 PtrList<volScalarField>& aLambda
329) const
330{
331 a = Zero;
332
333 for (label j=0; j<nBands_; j++)
334 {
335 aLambda[j].primitiveFieldRef() = this->a(j);
336
337 a.primitiveFieldRef() +=
338 aLambda[j].primitiveField()
339 *(iBands_[j][1] - iBands_[j][0])
340 /totalWaveLength_;
341 }
342
343}
344
345
346// ************************************************************************* //
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())
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
@ 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
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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
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.
wideBandAbsorptionEmission radiation absorption and emission coefficients for continuous phase.
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
FixedList< FixedList< absorptionCoeffs, nSpecies_ >, maxBands_ > coeffs_
Absorption coefficients.
tmp< volScalarField > aCont(const label bandi=0) const
Absorption coefficient for continuous phase.
wideBandAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
static const int nSpecies_
Maximum number of species considered for absorptivity.
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: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")
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.
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).
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
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.