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-2016 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  << "Model requires a multi-component thermo package"
71  << abort(FatalError);
72  }
73 
74 
75  label nFunc = 0;
76  const dictionary& functionDicts = dict.subDict(typeName + "Coeffs");
77 
78  forAllConstIter(dictionary, functionDicts, iter)
79  {
80  // safety:
81  if (!iter().isDict())
82  {
83  continue;
84  }
85  const word& key = iter().keyword();
86  speciesNames_.insert(key, nFunc);
87  const dictionary& dict = iter().dict();
88  coeffs_[nFunc].initialise(dict);
89  nFunc++;
90  }
91 
92  if (coeffsDict_.found("lookUpTableFileName"))
93  {
94  const word name = coeffsDict_.lookup("lookUpTableFileName");
95  if (name != "none")
96  {
97  lookUpTablePtr_.set
98  (
100  (
101  fileName(coeffsDict_.lookup("lookUpTableFileName")),
102  mesh.time().constant(),
103  mesh
104  )
105  );
106 
107  if (!mesh.foundObject<volScalarField>("ft"))
108  {
110  << "specie ft is not present to use with "
111  << "lookUpTableFileName " << nl
112  << exit(FatalError);
113  }
114  }
115  }
116 
117  // Check that all the species on the dictionary are present in the
118  // look-up table and save the corresponding indices of the look-up table
119 
120  label j = 0;
121  forAllConstIter(HashTable<label>, speciesNames_, iter)
122  {
123  if (!lookUpTablePtr_.empty())
124  {
125  if (lookUpTablePtr_().found(iter.key()))
126  {
127  label index = lookUpTablePtr_().findFieldIndex(iter.key());
128 
129  Info<< "specie: " << iter.key() << " found on look-up table "
130  << " with index: " << index << endl;
131 
132  specieIndex_[iter()] = index;
133  }
134  else if (mesh.foundObject<volScalarField>(iter.key()))
135  {
136  volScalarField& Y =
137  const_cast<volScalarField&>
138  (
139  mesh.lookupObject<volScalarField>(iter.key())
140  );
141  Yj_.set(j, &Y);
142  specieIndex_[iter()] = 0;
143  j++;
144  Info<< "specie: " << iter.key() << " is being solved" << endl;
145  }
146  else
147  {
149  << "specie: " << iter.key()
150  << " is neither in look-up table: "
151  << lookUpTablePtr_().tableName()
152  << " nor is being solved" << nl
153  << exit(FatalError);
154  }
155  }
156  else if (mesh.foundObject<volScalarField>(iter.key()))
157  {
158  volScalarField& Y =
159  const_cast<volScalarField&>
160  (
161  mesh.lookupObject<volScalarField>(iter.key())
162  );
163 
164  Yj_.set(j, &Y);
165  specieIndex_[iter()] = 0;
166  j++;
167  }
168  else
169  {
171  << " there is not lookup table and the specie" << nl
172  << iter.key() << nl
173  << " is not found " << nl
174  << exit(FatalError);
175 
176  }
177  }
178 }
179 
180 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
181 
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
190 {
191  const basicSpecieMixture& mixture =
192  dynamic_cast<const basicSpecieMixture&>(thermo_);
193 
194  const volScalarField& T = thermo_.T();
195  const volScalarField& p = thermo_.p();
196 
197 
199  (
200  new volScalarField
201  (
202  IOobject
203  (
204  "aCont" + name(bandI),
205  mesh().time().timeName(),
206  mesh(),
209  ),
210  mesh(),
212  extrapolatedCalculatedFvPatchVectorField::typeName
213  )
214  );
215 
216  scalarField& a = ta.ref().primitiveFieldRef();
217 
218  forAll(a, celli)
219  {
220  forAllConstIter(HashTable<label>, speciesNames_, iter)
221  {
222  label n = iter();
223  scalar Xipi = 0.0;
224  if (specieIndex_[n] != 0)
225  {
226  //Specie found in the lookUpTable.
227  const volScalarField& ft =
228  mesh_.lookupObject<volScalarField>("ft");
229 
230  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
231  //moles x pressure [atm]
232  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
233  }
234  else
235  {
236  scalar invWt = 0.0;
237  forAll(mixture.Y(), s)
238  {
239  invWt += mixture.Y(s)[celli]/mixture.W(s);
240  }
241 
242  label index = mixture.species()[iter.key()];
243  scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
244 
245  Xipi = Xk*paToAtm(p[celli]);
246  }
247 
248  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
249 
250  scalar Ti = T[celli];
251  // negative temperature exponents
252  if (coeffs_[n].invTemp())
253  {
254  Ti = 1.0/T[celli];
255  }
256  a[celli] +=
257  Xipi
258  *(
259  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
260  + b[0]
261  );
262  }
263  }
264  ta.ref().correctBoundaryConditions();
265  return ta;
266 }
267 
268 
271 {
272  return aCont(bandI);
273 }
274 
275 
278 {
280  (
281  new volScalarField
282  (
283  IOobject
284  (
285  "ECont" + name(bandI),
286  mesh_.time().timeName(),
287  mesh_,
290  ),
291  mesh_,
293  )
294  );
295 
296  if (mesh_.foundObject<volScalarField>("dQ"))
297  {
298  const volScalarField& dQ =
299  mesh_.lookupObject<volScalarField>("dQ");
300 
301  if (dQ.dimensions() == dimEnergy/dimTime)
302  {
303  E.ref().primitiveFieldRef() = EhrrCoeff_*dQ/mesh_.V();
304  }
305  else if (dQ.dimensions() == dimEnergy/dimTime/dimVolume)
306  {
307  E.ref().primitiveFieldRef() = EhrrCoeff_*dQ;
308  }
309  else
310  {
311  if (debug)
312  {
314  << "Incompatible dimensions for dQ field" << endl;
315  }
316  }
317  }
318  else
319  {
321  << "dQ field not found in mesh" << endl;
322  }
323 
324  return E;
325 }
326 
327 
328 // ************************************************************************* //
#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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool foundObject(const word &name) const
Is the named Type found?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Unit conversion functions.
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:633
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Model to supply absorption and emission coefficients for radiation modelling.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Macros for easy insertion into run-time selection tables.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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
dQ
Definition: YEEqn.H:14
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
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet & dimensions() const
Return dimensions.
static const char nl
Definition: Ostream.H:262
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
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
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.
PtrList< volScalarField > & Y
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
const speciesTable & species() const
Return the table of species.
label n
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:54
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243