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  speciesNames_(0),
62  specieIndex_(label(0)),
63  lookUpTablePtr_(),
64  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
65  Yj_(nSpecies_),
66  totalWaveLength_(0)
67 {
68  if (!isA<fluidMulticomponentThermo>(thermo_))
69  {
71  << "Model requires a multi-component thermo package"
72  << abort(FatalError);
73  }
74 
75  label nBand = 0;
77  {
78  if (!iter().isDict()) continue;
79 
80  const dictionary& dict = iter().dict();
81  dict.lookup("bandLimits") >> iBands_[nBand];
82  totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
83 
84  label nSpec = 0;
86  {
87  if (!iter().isDict()) continue;
88 
89  const word& key = iter().keyword();
90  if (nBand == 0)
91  {
92  speciesNames_.insert(key, nSpec);
93  }
94  else
95  {
96  if (!speciesNames_.found(key))
97  {
99  << "specie: " << key << " is not in all the bands"
100  << nl << exit(FatalError);
101  }
102  }
103  coeffs_[nBand][nSpec].initialise(dict.subDict(key));
104  nSpec++;
105  }
106  nBand++;
107  }
108  nBands_ = nBand;
109 
110  if (dict.found("lookUpTableFileName"))
111  {
112  const word name = dict.lookup("lookUpTableFileName");
113  if (name != "none")
114  {
115  lookUpTablePtr_.set
116  (
118  (
119  fileName(dict.lookup("lookUpTableFileName")),
120  mesh.time().constant(),
121  mesh
122  )
123  );
124 
125  if (!mesh.foundObject<volScalarField>("ft"))
126  {
128  << "specie ft is not present to use with "
129  << "lookUpTableFileName " << nl
130  << exit(FatalError);
131  }
132  }
133  }
134 
135  // Check that all the species on the dictionary are present in the
136  // look-up table and save the corresponding indices of the look-up table
137 
138  label j = 0;
140  {
141  if (!lookUpTablePtr_.empty())
142  {
143  if (lookUpTablePtr_().found(iter.key()))
144  {
145  const label index =
146  lookUpTablePtr_().findFieldIndex(iter.key());
147 
148  Info<< "specie: " << iter.key() << " found on look-up table "
149  << " with index: " << index << endl;
150 
151  specieIndex_[iter()] = index;
152  }
153  else if (mesh.foundObject<volScalarField>(iter.key()))
154  {
155  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
156  specieIndex_[iter()] = 0;
157  j++;
158  Info<< "specie: " << iter.key() << " is being solved" << endl;
159  }
160  else
161  {
163  << "specie: " << iter.key()
164  << " is neither in look-up table: "
165  << lookUpTablePtr_().tableName()
166  << " nor is being solved" << nl
167  << exit(FatalError);
168  }
169  }
170  else if (mesh.foundObject<volScalarField>(iter.key()))
171  {
172  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
173  specieIndex_[iter()] = 0;
174  j++;
175  }
176  else
177  {
179  << " there is no lookup table and the specie" << nl
180  << iter.key() << nl
181  << " is not found " << nl
182  << exit(FatalError);
183 
184  }
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
190 
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
199 (
200  const label bandi
201 ) const
202 {
203  const fluidMulticomponentThermo& mcThermo =
204  dynamic_cast<const fluidMulticomponentThermo&>(thermo_);
205 
206  const volScalarField& T = thermo_.T();
207  const volScalarField& p = thermo_.p();
208 
210  (
212  (
213  "a",
214  mesh(),
216  )
217  );
218 
219  scalarField& a = ta.ref().primitiveFieldRef();
220 
221  const unitConversion& unitAtm = units()["atm"];
222 
223  forAll(a, celli)
224  {
225  forAllConstIter(HashTable<label>, speciesNames_, iter)
226  {
227  const label n = iter();
228  scalar Xipi = 0;
229  if (specieIndex_[n] != 0)
230  {
231  const volScalarField& ft =
232  mesh_.lookupObject<volScalarField>("ft");
233 
234  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
235 
236  // moles*pressure [atm]
237  Xipi = unitAtm.toUser(Ynft[specieIndex_[n]]*p[celli]);
238  }
239  else
240  {
241  scalar invWt = 0;
242  forAll(mcThermo.Y(), s)
243  {
244  invWt += mcThermo.Y(s)[celli]/mcThermo.WiValue(s);
245  }
246 
247  const label index = mcThermo.species()[iter.key()];
248 
249  const scalar Xk =
250  mcThermo.Y(index)[celli]/(mcThermo.WiValue(index)*invWt);
251 
252  Xipi = unitAtm.toUser(Xk*p[celli]);
253  }
254 
255  scalar Ti = T[celli];
256 
258  coeffs_[bandi][n].coeffs(T[celli]);
259 
260  if (coeffs_[bandi][n].invTemp())
261  {
262  Ti = 1.0/T[celli];
263  }
264 
265  a[celli]+=
266  Xipi
267  *(
268  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
269  + b[0]
270  );
271  }
272  }
273 
274  return ta;
275 }
276 
277 
280 (
281  const label bandi
282 ) const
283 {
284  return aCont(bandi);
285 }
286 
287 
290 (
291  const label bandi
292 ) const
293 {
294  return absorptionEmissionModel::ECont(bandi);
295 }
296 
297 
299 (
300  volScalarField& a,
301  PtrList<volScalarField>& aLambda
302 ) const
303 {
305 
306  for (label j=0; j<nBands_; j++)
307  {
308  aLambda[j].primitiveFieldRef() = this->a(j);
309 
310  a.primitiveFieldRef() +=
311  aLambda[j].primitiveField()
312  *(iBands_[j][1] - iBands_[j][0])
313  /totalWaveLength_;
314  }
315 
316 }
317 
318 
319 // ************************************************************************* //
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
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 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:157
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
Definition: wideBand.C:299
UPtrList< volScalarField > Yj_
Pointer list of species being solved involved in the absorption.
Definition: wideBand.H:163
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:199
HashTable< label > speciesNames_
Hash table with species names.
Definition: wideBand.H:145
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:154
FixedList< label, nSpecies_ > specieIndex_
Indices of species in the look-up table.
Definition: wideBand.H:148
tmp< volScalarField > ECont(const label bandi=0) const
Emission contribution for continuous phase.
Definition: wideBand.C:290
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
Definition: wideBand.C:280
FixedList< Vector2D< scalar >, maxBands_ > iBands_
Bands.
Definition: wideBand.H:151
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:197
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
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:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
const dimensionSet dimLength
const HashTable< unitConversion > & units()
Get the table of unit conversions.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p