wideBandAbsorptionEmission.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-2018 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 
28 #include "basicSpecieMixture.H"
29 #include "unitConversion.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace radiation
36  {
37  defineTypeNameAndDebug(wideBandAbsorptionEmission, 0);
38 
40  (
41  absorptionEmissionModel,
42  wideBandAbsorptionEmission,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& dict,
54  const fvMesh& mesh
55 )
56 :
57  absorptionEmissionModel(dict, mesh),
58  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
59  speciesNames_(0),
60  specieIndex_(label(0)),
61  lookUpTablePtr_(),
63  Yj_(nSpecies_),
64  totalWaveLength_(0)
65 {
66  label nBand = 0;
67  const dictionary& functionDicts = dict.optionalSubDict(typeName +"Coeffs");
68  forAllConstIter(dictionary, functionDicts, iter)
69  {
70  // safety:
71  if (!iter().isDict())
72  {
73  continue;
74  }
75 
76  const dictionary& dict = iter().dict();
77  dict.lookup("bandLimits") >> iBands_[nBand];
78  dict.lookup("EhrrCoeff") >> iEhrrCoeffs_[nBand];
79  totalWaveLength_ += iBands_[nBand][1] - iBands_[nBand][0];
80 
81  label nSpec = 0;
82  const dictionary& specDicts = dict.subDict("species");
83  forAllConstIter(dictionary, specDicts, iter)
84  {
85  const word& key = iter().keyword();
86  if (nBand == 0)
87  {
88  speciesNames_.insert(key, nSpec);
89  }
90  else
91  {
92  if (!speciesNames_.found(key))
93  {
95  << "specie: " << key << " is not in all the bands"
96  << nl << exit(FatalError);
97  }
98  }
99  coeffs_[nBand][nSpec].initialise(specDicts.subDict(key));
100  nSpec++;
101  }
102  nBand++;
103  }
104  nBands_ = nBand;
105 
106  if (coeffsDict_.found("lookUpTableFileName"))
107  {
108  const word name = coeffsDict_.lookup("lookUpTableFileName");
109  if (name != "none")
110  {
111  lookUpTablePtr_.set
112  (
114  (
115  fileName(coeffsDict_.lookup("lookUpTableFileName")),
116  mesh.time().constant(),
117  mesh
118  )
119  );
120 
121  if (!mesh.foundObject<volScalarField>("ft"))
122  {
124  << "specie ft is not present to use with "
125  << "lookUpTableFileName " << nl
126  << exit(FatalError);
127  }
128  }
129  }
130 
131  // Check that all the species on the dictionary are present in the
132  // look-up table and save the corresponding indices of the look-up table
133 
134  label j = 0;
135  forAllConstIter(HashTable<label>, speciesNames_, iter)
136  {
137  if (!lookUpTablePtr_.empty())
138  {
139  if (lookUpTablePtr_().found(iter.key()))
140  {
141  const label index =
142  lookUpTablePtr_().findFieldIndex(iter.key());
143 
144  Info<< "specie: " << iter.key() << " found on look-up table "
145  << " with index: " << index << endl;
146 
147  specieIndex_[iter()] = index;
148  }
149  else if (mesh.foundObject<volScalarField>(iter.key()))
150  {
151  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
152  specieIndex_[iter()] = 0;
153  j++;
154  Info<< "specie: " << iter.key() << " is being solved" << endl;
155  }
156  else
157  {
159  << "specie: " << iter.key()
160  << " is neither in look-up table: "
161  << lookUpTablePtr_().tableName()
162  << " nor is being solved" << nl
163  << exit(FatalError);
164  }
165  }
166  else if (mesh.foundObject<volScalarField>(iter.key()))
167  {
168  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
169  specieIndex_[iter()] = 0;
170  j++;
171  }
172  else
173  {
175  << " there is no lookup table and the specie" << nl
176  << iter.key() << nl
177  << " is not found " << nl
178  << exit(FatalError);
179 
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
186 
188 {}
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
195 {
196  const basicSpecieMixture& mixture =
197  dynamic_cast<const basicSpecieMixture&>(thermo_);
198 
199  const volScalarField& T = thermo_.T();
200  const volScalarField& p = thermo_.p();
201 
203  (
204  new volScalarField
205  (
206  IOobject
207  (
208  "a",
209  mesh().time().timeName(),
210  mesh(),
213  ),
214  mesh(),
216  )
217  );
218 
219  scalarField& a = ta.ref().primitiveFieldRef();
220 
221  forAll(a, celli)
222  {
223  forAllConstIter(HashTable<label>, speciesNames_, iter)
224  {
225  const label n = iter();
226  scalar Xipi = 0;
227  if (specieIndex_[n] != 0)
228  {
229  const volScalarField& ft =
230  mesh_.lookupObject<volScalarField>("ft");
231 
232  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
233 
234  // moles*pressure [atm]
235  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
236  }
237  else
238  {
239  scalar invWt = 0;
240  forAll(mixture.Y(), s)
241  {
242  invWt += mixture.Y(s)[celli]/mixture.W(s);
243  }
244 
245  const label index = mixture.species()[iter.key()];
246 
247  const scalar Xk =
248  mixture.Y(index)[celli]/(mixture.W(index)*invWt);
249 
250  Xipi = Xk*paToAtm(p[celli]);
251  }
252 
253  scalar Ti = T[celli];
254 
256  coeffs_[bandi][n].coeffs(T[celli]);
257 
258  if (coeffs_[bandi][n].invTemp())
259  {
260  Ti = 1.0/T[celli];
261  }
262 
263  a[celli]+=
264  Xipi
265  *(
266  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
267  + b[0]
268  );
269  }
270  }
271 
272  return ta;
273 }
274 
275 
278 {
279  return aCont(bandi);
280 }
281 
282 
285 {
287  (
288  new volScalarField
289  (
290  IOobject
291  (
292  "E",
293  mesh().time().timeName(),
294  mesh(),
297  ),
298  mesh(),
300  )
301  );
302 
303  if (mesh().foundObject<volScalarField>("Qdot"))
304  {
305  const volScalarField& Qdot =
306  mesh().lookupObject<volScalarField>("Qdot");
307 
308  if (Qdot.dimensions() == dimEnergy/dimTime)
309  {
310  E.ref().primitiveFieldRef() =
311  iEhrrCoeffs_[bandi]
312  *Qdot.primitiveField()
313  *(iBands_[bandi][1] - iBands_[bandi][0])
314  /totalWaveLength_
315  /mesh_.V();
316  }
317  else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
318  {
319  E.ref().primitiveFieldRef() =
320  iEhrrCoeffs_[bandi]
321  *Qdot.primitiveField()
322  *(iBands_[bandi][1] - iBands_[bandi][0])
323  /totalWaveLength_;
324  }
325  else
326  {
328  << "Incompatible dimensions for Qdot field" << endl;
329  }
330  }
331 
332  return E;
333 }
334 
335 
337 (
338  volScalarField& a,
339  PtrList<volScalarField>& aLambda
340 ) const
341 {
342  a = dimensionedScalar("zero", dimless/dimLength, 0);
343 
344  for (label j=0; j<nBands_; j++)
345  {
346  aLambda[j].primitiveFieldRef() = this->a(j);
347 
348  a.primitiveFieldRef() +=
349  aLambda[j].primitiveField()
350  *(iBands_[j][1] - iBands_[j][0])
351  /totalWaveLength_;
352  }
353 
354 }
355 
356 
357 // ************************************************************************* //
tmp< volScalarField > aCont(const label bandi=0) const
Absorption coefficient for continuous phase.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void correct(volScalarField &a, PtrList< volScalarField > &aLambda) const
Correct rays.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Unit conversion functions.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const speciesTable & species() const
Return the table of species.
bool foundObject(const word &name) const
Is the named Type found?
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Model to supply absorption and emission coefficients for radiation modelling.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
scalar Qdot
Definition: solveChemistry.H:2
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:49
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
word timeName
Definition: getTimeIndex.H:3
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:265
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
label n
wideBandAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< volScalarField > ECont(const label bandi=0) const
Emission contribution for continuous phase.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
tmp< volScalarField > eCont(const label bandi=0) const
Emission coefficient for continuous phase.