wideBand.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wideBand.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace radiationModels
35 {
36 namespace absorptionEmissionModels
37 {
39 
41  (
43  wideBand,
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict,
56  const fvMesh& mesh,
57  const word& modelName
58 )
59 :
61  coeffsDict_(dict.subDict(modelName + "Coeffs")),
62  speciesNames_(0),
63  specieIndex_(label(0)),
64  lookUpTablePtr_(),
65  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
66  Yj_(nSpecies_),
67  totalWaveLength_(0)
68 {
69  if (!isA<fluidMulticomponentThermo>(thermo_))
70  {
72  << "Model requires a multi-component thermo package"
73  << abort(FatalError);
74  }
75 
76  label nBand = 0;
78  {
79  if (!iter().isDict()) continue;
80 
81  const dictionary& dict = iter().dict();
82  dict.lookup("bandLimits") >> iBands_[nBand];
83  totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
84 
85  label nSpec = 0;
87  {
88  if (!iter().isDict()) continue;
89 
90  const word& key = iter().keyword();
91  if (nBand == 0)
92  {
93  speciesNames_.insert(key, nSpec);
94  }
95  else
96  {
97  if (!speciesNames_.found(key))
98  {
100  << "specie: " << key << " is not in all the bands"
101  << nl << exit(FatalError);
102  }
103  }
104  coeffs_[nBand][nSpec].initialise(dict.subDict(key));
105  nSpec++;
106  }
107  nBand++;
108  }
109  nBands_ = nBand;
110 
111  if (coeffsDict_.found("lookUpTableFileName"))
112  {
113  const word name = coeffsDict_.lookup("lookUpTableFileName");
114  if (name != "none")
115  {
116  lookUpTablePtr_.set
117  (
119  (
120  fileName(coeffsDict_.lookup("lookUpTableFileName")),
121  mesh.time().constant(),
122  mesh
123  )
124  );
125 
126  if (!mesh.foundObject<volScalarField>("ft"))
127  {
129  << "specie ft is not present to use with "
130  << "lookUpTableFileName " << nl
131  << exit(FatalError);
132  }
133  }
134  }
135 
136  // Check that all the species on the dictionary are present in the
137  // look-up table and save the corresponding indices of the look-up table
138 
139  label j = 0;
141  {
142  if (!lookUpTablePtr_.empty())
143  {
144  if (lookUpTablePtr_().found(iter.key()))
145  {
146  const label index =
147  lookUpTablePtr_().findFieldIndex(iter.key());
148 
149  Info<< "specie: " << iter.key() << " found on look-up table "
150  << " with index: " << index << endl;
151 
152  specieIndex_[iter()] = index;
153  }
154  else if (mesh.foundObject<volScalarField>(iter.key()))
155  {
156  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
157  specieIndex_[iter()] = 0;
158  j++;
159  Info<< "specie: " << iter.key() << " is being solved" << endl;
160  }
161  else
162  {
164  << "specie: " << iter.key()
165  << " is neither in look-up table: "
166  << lookUpTablePtr_().tableName()
167  << " nor is being solved" << nl
168  << exit(FatalError);
169  }
170  }
171  else if (mesh.foundObject<volScalarField>(iter.key()))
172  {
173  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
174  specieIndex_[iter()] = 0;
175  j++;
176  }
177  else
178  {
180  << " there is no lookup table and the specie" << nl
181  << iter.key() << nl
182  << " is not found " << nl
183  << exit(FatalError);
184 
185  }
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 
193 {}
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
200 (
201  const label bandi
202 ) const
203 {
204  const fluidMulticomponentThermo& mcThermo =
205  dynamic_cast<const fluidMulticomponentThermo&>(thermo_);
206 
207  const volScalarField& T = thermo_.T();
208  const volScalarField& p = thermo_.p();
209 
211  (
213  (
214  "a",
215  mesh(),
217  )
218  );
219 
220  scalarField& a = ta.ref().primitiveFieldRef();
221 
222  const unitConversion& unitAtm = units()["atm"];
223 
224  forAll(a, celli)
225  {
226  forAllConstIter(HashTable<label>, speciesNames_, iter)
227  {
228  const label n = iter();
229  scalar Xipi = 0;
230  if (specieIndex_[n] != 0)
231  {
232  const volScalarField& ft =
233  mesh_.lookupObject<volScalarField>("ft");
234 
235  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
236 
237  // moles*pressure [atm]
238  Xipi = unitAtm.toUser(Ynft[specieIndex_[n]]*p[celli]);
239  }
240  else
241  {
242  scalar invWt = 0;
243  forAll(mcThermo.Y(), s)
244  {
245  invWt += mcThermo.Y(s)[celli]/mcThermo.WiValue(s);
246  }
247 
248  const label index = mcThermo.species()[iter.key()];
249 
250  const scalar Xk =
251  mcThermo.Y(index)[celli]/(mcThermo.WiValue(index)*invWt);
252 
253  Xipi = unitAtm.toUser(Xk*p[celli]);
254  }
255 
256  scalar Ti = T[celli];
257 
259  coeffs_[bandi][n].coeffs(T[celli]);
260 
261  if (coeffs_[bandi][n].invTemp())
262  {
263  Ti = 1.0/T[celli];
264  }
265 
266  a[celli]+=
267  Xipi
268  *(
269  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
270  + b[0]
271  );
272  }
273  }
274 
275  return ta;
276 }
277 
278 
281 (
282  const label bandi
283 ) const
284 {
285  return aCont(bandi);
286 }
287 
288 
291 (
292  const label bandi
293 ) const
294 {
295  return absorptionEmissionModel::ECont(bandi);
296 }
297 
298 
300 (
301  volScalarField& a,
302  PtrList<volScalarField>& aLambda
303 ) const
304 {
306 
307  for (label j=0; j<nBands_; j++)
308  {
309  aLambda[j].primitiveFieldRef() = this->a(j);
310 
311  a.primitiveFieldRef() +=
312  aLambda[j].primitiveField()
313  *(iBands_[j][1] - iBands_[j][0])
314  /totalWaveLength_;
315  }
316 
317 }
318 
319 
320 // ************************************************************************* //
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
A class for handling file names.
Definition: fileName.H:82
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
virtual const speciesTable & species() const =0
The table of species.
virtual scalar WiValue(const label speciei) const =0
Molecular weight [kg/kmol].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
bool foundObject(const word &name) const
Is the named Type in registry.
A base class for physical properties.
Model to supply absorption and emission coefficients for radiation modelling.
const dictionary & dict() const
Reference to the dictionary.
const fvMesh & mesh() const
Reference to the mesh.
virtual tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
const fluidThermo & thermo_
Thermo package.
Definition: wideBand.H:160
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
Definition: wideBand.C:300
UPtrList< volScalarField > Yj_
Pointer list of species being solved involved in the absorption.
Definition: wideBand.H:166
FixedList< FixedList< absorptionCoeffs, nSpecies_ >, maxBands_ > coeffs_
Absorption coefficients.
Definition: wideBand.H:137
tmp< volScalarField > aCont(const label bandi=0) const
Absorption coefficient for continuous phase.
Definition: wideBand.C:200
HashTable< label > speciesNames_
Hash table with species names.
Definition: wideBand.H:148
wideBand(const dictionary &dict, const fvMesh &mesh, const word &modelName=typeName)
Construct from components.
Definition: wideBand.C:54
autoPtr< interpolationLookUpTable > lookUpTablePtr_
Look-up table of species related to ft.
Definition: wideBand.H:157
FixedList< label, nSpecies_ > specieIndex_
Indices of species in the look-up table.
Definition: wideBand.H:151
tmp< volScalarField > ECont(const label bandi=0) const
Emission contribution for continuous phase.
Definition: wideBand.C:291
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
Definition: wideBand.C:281
FixedList< Vector2D< scalar >, maxBands_ > iBands_
Bands.
Definition: wideBand.H:154
dictionary coeffsDict_
Absorption model dictionary.
Definition: wideBand.H:145
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
T toUser(const T &) const
Convert a value to user units.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:25
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
const HashTable< unitConversion > & units()
Get the table of unit conversions.
error FatalError
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p