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-2017 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.optionalSubDict(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.optionalSubDict(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  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
137  specieIndex_[iter()] = 0;
138  j++;
139  Info<< "specie: " << iter.key() << " is being solved" << endl;
140  }
141  else
142  {
144  << "specie: " << iter.key()
145  << " is neither in look-up table: "
146  << lookUpTablePtr_().tableName()
147  << " nor is being solved" << nl
148  << exit(FatalError);
149  }
150  }
151  else if (mesh.foundObject<volScalarField>(iter.key()))
152  {
153  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
154  specieIndex_[iter()] = 0;
155  j++;
156  }
157  else
158  {
160  << " there is not lookup table and the specie" << nl
161  << iter.key() << nl
162  << " is not found " << nl
163  << exit(FatalError);
164 
165  }
166  }
167 }
168 
169 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
170 
172 {}
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
179 {
180  const basicSpecieMixture& mixture =
181  dynamic_cast<const basicSpecieMixture&>(thermo_);
182 
183  const volScalarField& T = thermo_.T();
184  const volScalarField& p = thermo_.p();
185 
186 
188  (
189  new volScalarField
190  (
191  IOobject
192  (
193  "aCont" + name(bandI),
194  mesh().time().timeName(),
195  mesh(),
198  ),
199  mesh(),
201  extrapolatedCalculatedFvPatchVectorField::typeName
202  )
203  );
204 
205  scalarField& a = ta.ref().primitiveFieldRef();
206 
207  forAll(a, celli)
208  {
209  forAllConstIter(HashTable<label>, speciesNames_, iter)
210  {
211  label n = iter();
212  scalar Xipi = 0.0;
213  if (specieIndex_[n] != 0)
214  {
215  //Specie found in the lookUpTable.
216  const volScalarField& ft =
217  mesh_.lookupObject<volScalarField>("ft");
218 
219  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
220  //moles x pressure [atm]
221  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
222  }
223  else
224  {
225  scalar invWt = 0.0;
226  forAll(mixture.Y(), s)
227  {
228  invWt += mixture.Y(s)[celli]/mixture.W(s);
229  }
230 
231  label index = mixture.species()[iter.key()];
232  scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
233 
234  Xipi = Xk*paToAtm(p[celli]);
235  }
236 
237  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
238 
239  scalar Ti = T[celli];
240  // negative temperature exponents
241  if (coeffs_[n].invTemp())
242  {
243  Ti = 1.0/T[celli];
244  }
245  a[celli] +=
246  Xipi
247  *(
248  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
249  + b[0]
250  );
251  }
252  }
253  ta.ref().correctBoundaryConditions();
254  return ta;
255 }
256 
257 
260 {
261  return aCont(bandI);
262 }
263 
264 
267 {
269  (
270  new volScalarField
271  (
272  IOobject
273  (
274  "ECont" + name(bandI),
275  mesh_.time().timeName(),
276  mesh_,
279  ),
280  mesh_,
282  )
283  );
284 
285  if (mesh_.foundObject<volScalarField>("Qdot"))
286  {
287  const volScalarField& Qdot =
288  mesh_.lookupObject<volScalarField>("Qdot");
289 
290  if (Qdot.dimensions() == dimEnergy/dimTime)
291  {
292  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot/mesh_.V();
293  }
294  else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
295  {
296  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot;
297  }
298  else
299  {
300  if (debug)
301  {
303  << "Incompatible dimensions for Qdot field" << endl;
304  }
305  }
306  }
307  else
308  {
310  << "Qdot field not found in mesh" << endl;
311  }
312 
313  return E;
314 }
315 
316 
317 // ************************************************************************* //
#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
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
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
#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:253
const speciesTable & species() const
Return the table of species.
bool foundObject(const word &name) const
Is the named Type found?
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.
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
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 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
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
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
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
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.
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
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
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.