greyMean.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-2024 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 
26 #include "greyMean.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace radiationModels
36 {
37 namespace absorptionEmissionModels
38 {
40 
42  (
44  greyMean,
46  );
47 }
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const dictionary& dict,
57  const fvMesh& mesh,
58  const word& modelName
59 )
60 :
62  coeffsDict_(dict.subDict(modelName + "Coeffs")),
63  speciesNames_(0),
64  specieIndex_(label(0)),
65  lookUpTablePtr_(),
66  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
67  Yj_(nSpecies_)
68 {
69  if (!isA<fluidMulticomponentThermo>(thermo_))
70  {
72  << "Model requires a multi-component thermo package"
73  << abort(FatalError);
74  }
75 
76  label nFunc = 0;
78  {
79  // safety:
80  if (!iter().isDict())
81  {
82  continue;
83  }
84  const word& key = iter().keyword();
85  speciesNames_.insert(key, nFunc);
86  const dictionary& dict = iter().dict();
87  coeffs_[nFunc].initialise(dict);
88  nFunc++;
89  }
90 
91  if (coeffsDict_.found("lookUpTableFileName"))
92  {
93  const word name = coeffsDict_.lookup("lookUpTableFileName");
94  if (name != "none")
95  {
96  lookUpTablePtr_.set
97  (
99  (
100  fileName(coeffsDict_.lookup("lookUpTableFileName")),
101  mesh.time().constant(),
102  mesh
103  )
104  );
105 
106  if (!mesh.foundObject<volScalarField>("ft"))
107  {
109  << "specie ft is not present to use with "
110  << "lookUpTableFileName " << nl
111  << exit(FatalError);
112  }
113  }
114  }
115 
116  // Check that all the species on the dictionary are present in the
117  // look-up table and save the corresponding indices of the look-up table
118 
119  label j = 0;
121  {
122  if (!lookUpTablePtr_.empty())
123  {
124  if (lookUpTablePtr_().found(iter.key()))
125  {
126  label index = lookUpTablePtr_().findFieldIndex(iter.key());
127 
128  Info<< "specie: " << iter.key() << " found on look-up table "
129  << " with index: " << index << endl;
130 
131  specieIndex_[iter()] = index;
132  }
133  else if (mesh.foundObject<volScalarField>(iter.key()))
134  {
135  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
136  specieIndex_[iter()] = 0;
137  j++;
138  Info<< "specie: " << iter.key() << " is being solved" << endl;
139  }
140  else
141  {
143  << "specie: " << iter.key()
144  << " is neither in look-up table: "
145  << lookUpTablePtr_().tableName()
146  << " nor is being solved" << nl
147  << exit(FatalError);
148  }
149  }
150  else if (mesh.foundObject<volScalarField>(iter.key()))
151  {
152  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
153  specieIndex_[iter()] = 0;
154  j++;
155  }
156  else
157  {
159  << " there is not lookup table and the specie" << nl
160  << iter.key() << nl
161  << " is not found " << nl
162  << exit(FatalError);
163 
164  }
165  }
166 }
167 
168 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
169 
171 {}
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
178 (
179  const label bandI
180 ) const
181 {
182  const fluidMulticomponentThermo& mcThermo =
183  dynamic_cast<const fluidMulticomponentThermo&>(thermo_);
184 
185  const volScalarField& T = thermo_.T();
186  const volScalarField& p = thermo_.p();
187 
188 
190  (
192  (
193  "aCont" + name(bandI),
194  mesh(),
196  extrapolatedCalculatedFvPatchVectorField::typeName
197  )
198  );
199 
200  scalarField& a = ta.ref().primitiveFieldRef();
201 
202  const unitConversion& unitAtm = units()["atm"];
203 
204  forAll(a, celli)
205  {
206  forAllConstIter(HashTable<label>, speciesNames_, iter)
207  {
208  label n = iter();
209  scalar Xipi = 0.0;
210  if (specieIndex_[n] != 0)
211  {
212  // Specie found in the lookUpTable.
213  const volScalarField& ft =
214  mesh_.lookupObject<volScalarField>("ft");
215 
216  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
217 
218  // moles x pressure [atm]
219  Xipi = unitAtm.toUser(Ynft[specieIndex_[n]]*p[celli]);
220  }
221  else
222  {
223  scalar invWt = 0.0;
224  forAll(mcThermo.Y(), s)
225  {
226  invWt += mcThermo.Y(s)[celli]/mcThermo.WiValue(s);
227  }
228 
229  const label index = mcThermo.species()[iter.key()];
230 
231  const scalar Xk =
232  mcThermo.Y(index)[celli]/(mcThermo.WiValue(index)*invWt);
233 
234  Xipi = unitAtm.toUser(Xk*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  const label bandI
262 ) const
263 {
264  return aCont(bandI);
265 }
266 
267 
270 (
271  const label bandI
272 ) const
273 {
274  return absorptionEmissionModel::ECont(bandI);
275 }
276 
277 
278 // ************************************************************************* //
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
A class for handling file names.
Definition: fileName.H:82
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
virtual const speciesTable & species() const =0
The table of species.
virtual scalar WiValue(const label speciei) const =0
Molecular weight [kg/kmol].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
bool foundObject(const word &name) const
Is the named Type in registry.
A base class for physical properties.
Model to supply absorption and emission coefficients for radiation modelling.
const dictionary & dict() const
Reference to the dictionary.
const fvMesh & mesh() const
Reference to the mesh.
virtual tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
const fluidThermo & thermo_
Thermo package.
Definition: greyMean.H:142
UPtrList< volScalarField > Yj_
Pointer list of species in the registry involved in the absorption.
Definition: greyMean.H:145
HashTable< label > speciesNames_
Hash table of species names.
Definition: greyMean.H:133
autoPtr< interpolationLookUpTable > lookUpTablePtr_
Look-up table of species related to ft.
Definition: greyMean.H:139
FixedList< label, nSpecies_ > specieIndex_
Indices of species in the look-up table.
Definition: greyMean.H:136
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMean.C:260
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMean.C:178
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
Definition: greyMean.C:270
dictionary coeffsDict_
Absorption model dictionary.
Definition: greyMean.H:130
greyMean(const dictionary &dict, const fvMesh &mesh, const word &modelName=typeName)
Construct from components.
Definition: greyMean.C:55
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
T toUser(const T &) const
Convert a value to user units.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:25
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
const HashTable< unitConversion > & units()
Get the table of unit conversions.
error FatalError
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p