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