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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "unitConversion.H"
33 #include "basicSpecieMixture.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace radiation
40  {
41  defineTypeNameAndDebug(greyMeanAbsorptionEmission, 0);
42 
44  (
45  absorptionEmissionModel,
46  greyMeanAbsorptionEmission,
47  dictionary
48  );
49  }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const dictionary& dict,
58  const fvMesh& mesh
59 )
60 :
62  coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
63  speciesNames_(0),
64  specieIndex_(Zero),
65  lookUpTablePtr_(),
66  thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
67  EhrrCoeff_(coeffsDict_.get<scalar>("EhrrCoeff")),
68  Yj_(nSpecies_)
69 {
70  if (!isA<basicSpecieMixture>(thermo_))
71  {
73  << "Model requires a multi-component thermo package"
74  << abort(FatalError);
75  }
76 
77 
78  label nFunc = 0;
79  const dictionary& functionDicts = dict.optionalSubDict(typeName + "Coeffs");
80 
81  for (const entry& dEntry : functionDicts)
82  {
83  if (!dEntry.isDict()) // safety
84  {
85  continue;
86  }
87 
88  const word& key = dEntry.keyword();
89  const dictionary& dict = dEntry.dict();
90 
91  speciesNames_.insert(key, nFunc);
92 
93  coeffs_[nFunc].initialise(dict);
94  nFunc++;
95  }
96 
97  if
98  (
99  coeffsDict_.found("lookUpTableFileName")
100  && "none" != coeffsDict_.get<word>("lookUpTableFileName")
101  )
102  {
103  lookUpTablePtr_.reset
104  (
105  new interpolationLookUpTable<scalar>
106  (
107  coeffsDict_.get<fileName>("lookUpTableFileName"),
108  mesh.time().constant(),
109  mesh
110  )
111  );
112 
113  if (!mesh.foundObject<volScalarField>("ft"))
114  {
116  << "specie ft is not present to use with "
117  << "lookUpTableFileName " << nl
118  << exit(FatalError);
119  }
120  }
121 
122  // Check that all the species on the dictionary are present in the
123  // look-up table and save the corresponding indices of the look-up table
124 
125  label j = 0;
126  forAllConstIters(speciesNames_, iter)
127  {
128  const word& specieName = iter.key();
129  const label index = iter.val();
130 
131  volScalarField* fldPtr = mesh.getObjectPtr<volScalarField>(specieName);
132 
133  if (lookUpTablePtr_)
134  {
135  if (lookUpTablePtr_().found(specieName))
136  {
137  const label fieldIndex =
138  lookUpTablePtr_().findFieldIndex(specieName);
139 
140  Info<< "specie: " << specieName << " found on look-up table "
141  << " with index: " << fieldIndex << endl;
142 
143  specieIndex_[index] = fieldIndex;
144  }
145  else if (fldPtr)
146  {
147  Yj_.set(j, fldPtr);
148  specieIndex_[index] = 0;
149  j++;
150  Info<< "specie: " << iter.key() << " is being solved" << endl;
151  }
152  else
153  {
155  << "specie: " << specieName
156  << " is neither in look-up table: "
157  << lookUpTablePtr_().tableName()
158  << " nor is being solved" << nl
159  << exit(FatalError);
160  }
161  }
162  else if (fldPtr)
163  {
164  Yj_.set(j, fldPtr);
165  specieIndex_[index] = 0;
166  j++;
167  }
168  else
169  {
171  << "There is no lookup table and the specie" << nl
172  << specieName << nl
173  << " is not found " << nl
175  }
176  }
177 }
178 
179 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
180 
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
189 {
190  const basicSpecieMixture& mixture =
191  dynamic_cast<const basicSpecieMixture&>(thermo_);
192 
193  const volScalarField& T = thermo_.T();
194  const volScalarField& p = thermo_.p();
195 
196 
198  (
199  new volScalarField
200  (
201  IOobject
202  (
203  "aCont" + name(bandI),
204  mesh().time().timeName(),
205  mesh(),
208  ),
209  mesh(),
212  )
213  );
214 
215  scalarField& a = ta.ref().primitiveFieldRef();
216 
217  forAll(a, celli)
218  {
219  forAllConstIters(speciesNames_, iter)
220  {
221  label n = iter();
222  scalar Xipi = 0.0;
223  if (specieIndex_[n] != 0)
224  {
225  //Specie found in the lookUpTable.
226  const volScalarField& ft =
227  mesh_.lookupObject<volScalarField>("ft");
228 
229  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
230  //moles x pressure [atm]
231  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
232  }
233  else
234  {
235  scalar invWt = 0.0;
236  forAll(mixture.Y(), s)
237  {
238  invWt += mixture.Y(s)[celli]/mixture.W(s);
239  }
240 
241  const label index = mixture.species().find(iter.key());
242  scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
243 
244  Xipi = Xk*paToAtm(p[celli]);
245  }
246 
247  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
248 
249  scalar Ti = T[celli];
250  // negative temperature exponents
251  if (coeffs_[n].invTemp())
252  {
253  Ti = 1.0/T[celli];
254  }
255  a[celli] +=
256  Xipi
257  *(
258  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
259  + b[0]
260  );
261  }
262  }
263  ta.ref().correctBoundaryConditions();
264  return ta;
265 }
266 
267 
270 {
271  return aCont(bandI);
272 }
273 
274 
277 {
279  (
280  new volScalarField
281  (
282  IOobject
283  (
284  "ECont" + name(bandI),
285  mesh_.time().timeName(),
286  mesh_,
289  ),
290  mesh_,
292  )
293  );
294 
295  const volScalarField* QdotPtr = mesh_.findObject<volScalarField>("Qdot");
296 
297  if (QdotPtr)
298  {
299  const volScalarField& Qdot = *QdotPtr;
300 
301  if (Qdot.dimensions() == dimEnergy/dimTime)
302  {
303  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot/mesh_.V();
304  }
305  else if (Qdot.dimensions() == dimEnergy/dimTime/dimVolume)
306  {
307  E.ref().primitiveFieldRef() = EhrrCoeff_*Qdot;
308  }
309  else
310  {
311  if (debug)
312  {
314  << "Incompatible dimensions for Qdot field" << endl;
315  }
316  }
317  }
318  else
319  {
321  << "Qdot field not found in mesh" << endl;
322  }
323 
324  return E;
325 }
326 
327 
328 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
void initialise(const dictionary &)
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:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
Unit conversion functions.
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
const word dictName("faMeshDefinition")
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const fvMesh & mesh() const
Reference to the mesh.
const dictionary & dict() const
Reference to the dictionary.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Model to supply absorption and emission coefficients for radiation modelling.
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
scalar Qdot
Definition: solveChemistry.H:2
Macros for easy insertion into run-time selection tables.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
word timeName
Definition: getTimeIndex.H:3
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:560
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
greyMeanAbsorptionEmission(const dictionary &dict, const fvMesh &mesh)
Construct from components.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
#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.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
FixedList< scalar, nCoeffs_ > coeffArray
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127