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-2021 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"
28 #include "basicSpecieMixture.H"
29 #include "unitConversion.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace radiationModels
36 {
37 namespace absorptionEmissionModels
38 {
40 
42  (
44  wideBand,
46  );
47 }
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const dictionary& dict,
57  const fvMesh& mesh,
58  const word& modelName
59 )
60 :
62  coeffsDict_(dict.subDict(modelName + "Coeffs")),
63  speciesNames_(0),
64  specieIndex_(label(0)),
65  lookUpTablePtr_(),
66  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
67  Yj_(nSpecies_),
68  totalWaveLength_(0)
69 {
70  label nBand = 0;
72  {
73  if (!iter().isDict()) continue;
74 
75  const dictionary& dict = iter().dict();
76  dict.lookup("bandLimits") >> iBands_[nBand];
77  totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
78 
79  label nSpec = 0;
81  {
82  if (!iter().isDict()) continue;
83 
84  const word& key = iter().keyword();
85  if (nBand == 0)
86  {
87  speciesNames_.insert(key, nSpec);
88  }
89  else
90  {
91  if (!speciesNames_.found(key))
92  {
94  << "specie: " << key << " is not in all the bands"
95  << nl << exit(FatalError);
96  }
97  }
98  coeffs_[nBand][nSpec].initialise(dict.subDict(key));
99  nSpec++;
100  }
101  nBand++;
102  }
103  nBands_ = nBand;
104 
105  if (coeffsDict_.found("lookUpTableFileName"))
106  {
107  const word name = coeffsDict_.lookup("lookUpTableFileName");
108  if (name != "none")
109  {
110  lookUpTablePtr_.set
111  (
113  (
114  fileName(coeffsDict_.lookup("lookUpTableFileName")),
115  mesh.time().constant(),
116  mesh
117  )
118  );
119 
120  if (!mesh.foundObject<volScalarField>("ft"))
121  {
123  << "specie ft is not present to use with "
124  << "lookUpTableFileName " << nl
125  << exit(FatalError);
126  }
127  }
128  }
129 
130  // Check that all the species on the dictionary are present in the
131  // look-up table and save the corresponding indices of the look-up table
132 
133  label j = 0;
135  {
136  if (!lookUpTablePtr_.empty())
137  {
138  if (lookUpTablePtr_().found(iter.key()))
139  {
140  const label index =
141  lookUpTablePtr_().findFieldIndex(iter.key());
142 
143  Info<< "specie: " << iter.key() << " found on look-up table "
144  << " with index: " << index << endl;
145 
146  specieIndex_[iter()] = index;
147  }
148  else if (mesh.foundObject<volScalarField>(iter.key()))
149  {
150  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
151  specieIndex_[iter()] = 0;
152  j++;
153  Info<< "specie: " << iter.key() << " is being solved" << endl;
154  }
155  else
156  {
158  << "specie: " << iter.key()
159  << " is neither in look-up table: "
160  << lookUpTablePtr_().tableName()
161  << " nor is being solved" << nl
162  << exit(FatalError);
163  }
164  }
165  else if (mesh.foundObject<volScalarField>(iter.key()))
166  {
167  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
168  specieIndex_[iter()] = 0;
169  j++;
170  }
171  else
172  {
174  << " there is no lookup table and the specie" << nl
175  << iter.key() << nl
176  << " is not found " << nl
177  << exit(FatalError);
178 
179  }
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
194 (
195  const label bandi
196 ) const
197 {
198  const basicSpecieMixture& mixture =
199  dynamic_cast<const basicSpecieMixture&>(thermo_);
200 
201  const volScalarField& T = thermo_.T();
202  const volScalarField& p = thermo_.p();
203 
205  (
207  (
208  "a",
209  mesh(),
211  )
212  );
213 
214  scalarField& a = ta.ref().primitiveFieldRef();
215 
216  forAll(a, celli)
217  {
218  forAllConstIter(HashTable<label>, speciesNames_, iter)
219  {
220  const label n = iter();
221  scalar Xipi = 0;
222  if (specieIndex_[n] != 0)
223  {
224  const volScalarField& ft =
225  mesh_.lookupObject<volScalarField>("ft");
226 
227  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
228 
229  // moles*pressure [atm]
230  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
231  }
232  else
233  {
234  scalar invWt = 0;
235  forAll(mixture.Y(), s)
236  {
237  invWt += mixture.Y(s)[celli]/mixture.Wi(s);
238  }
239 
240  const label index = mixture.species()[iter.key()];
241 
242  const scalar Xk =
243  mixture.Y(index)[celli]/(mixture.Wi(index)*invWt);
244 
245  Xipi = Xk*paToAtm(p[celli]);
246  }
247 
248  scalar Ti = T[celli];
249 
251  coeffs_[bandi][n].coeffs(T[celli]);
252 
253  if (coeffs_[bandi][n].invTemp())
254  {
255  Ti = 1.0/T[celli];
256  }
257 
258  a[celli]+=
259  Xipi
260  *(
261  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
262  + b[0]
263  );
264  }
265  }
266 
267  return ta;
268 }
269 
270 
273 (
274  const label bandi
275 ) const
276 {
277  return aCont(bandi);
278 }
279 
280 
283 (
284  const label bandi
285 ) const
286 {
287  return absorptionEmissionModel::ECont(bandi);
288 }
289 
290 
292 (
293  volScalarField& a,
294  PtrList<volScalarField>& aLambda
295 ) const
296 {
298 
299  for (label j=0; j<nBands_; j++)
300  {
301  aLambda[j].primitiveFieldRef() = this->a(j);
302 
303  a.primitiveFieldRef() +=
304  aLambda[j].primitiveField()
305  *(iBands_[j][1] - iBands_[j][0])
306  /totalWaveLength_;
307  }
308 
309 }
310 
311 
312 // ************************************************************************* //
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 internal field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
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:113
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
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
A class for handling file names.
Definition: fileName.H:82
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
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.
An abstract 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.
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
Definition: wideBand.C:292
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:194
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:55
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:283
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.
Definition: wideBand.C:273
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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:27
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
error FatalError
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p
Unit conversion functions.