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  speciesNames_(0),
63  specieIndex_(label(0)),
64  lookUpTablePtr_(),
65  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
66  Yj_(nSpecies_)
67 {
68  if (!isA<fluidMulticomponentThermo>(thermo_))
69  {
71  << "Model requires a multi-component thermo package"
72  << abort(FatalError);
73  }
74 
75  label nFunc = 0;
77  {
78  // safety:
79  if (!iter().isDict())
80  {
81  continue;
82  }
83  const word& key = iter().keyword();
84  speciesNames_.insert(key, nFunc);
85  const dictionary& dict = iter().dict();
86  coeffs_[nFunc].initialise(dict);
87  nFunc++;
88  }
89 
90  if (dict.found("lookUpTableFileName"))
91  {
92  const word name = dict.lookup("lookUpTableFileName");
93  if (name != "none")
94  {
95  lookUpTablePtr_.set
96  (
98  (
99  fileName(dict.lookup("lookUpTableFileName")),
100  mesh.time().constant(),
101  mesh
102  )
103  );
104 
105  if (!mesh.foundObject<volScalarField>("ft"))
106  {
108  << "specie ft is not present to use with "
109  << "lookUpTableFileName " << nl
110  << exit(FatalError);
111  }
112  }
113  }
114 
115  // Check that all the species on the dictionary are present in the
116  // look-up table and save the corresponding indices of the look-up table
117 
118  label j = 0;
120  {
121  if (!lookUpTablePtr_.empty())
122  {
123  if (lookUpTablePtr_().found(iter.key()))
124  {
125  label index = lookUpTablePtr_().findFieldIndex(iter.key());
126 
127  Info<< "specie: " << iter.key() << " found on look-up table "
128  << " with index: " << index << endl;
129 
130  specieIndex_[iter()] = index;
131  }
132  else if (mesh.foundObject<volScalarField>(iter.key()))
133  {
134  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
135  specieIndex_[iter()] = 0;
136  j++;
137  Info<< "specie: " << iter.key() << " is being solved" << endl;
138  }
139  else
140  {
142  << "specie: " << iter.key()
143  << " is neither in look-up table: "
144  << lookUpTablePtr_().tableName()
145  << " nor is being solved" << nl
146  << exit(FatalError);
147  }
148  }
149  else if (mesh.foundObject<volScalarField>(iter.key()))
150  {
151  Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
152  specieIndex_[iter()] = 0;
153  j++;
154  }
155  else
156  {
158  << " there is not lookup table and the specie" << nl
159  << iter.key() << nl
160  << " is not found " << nl
161  << exit(FatalError);
162 
163  }
164  }
165 }
166 
167 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
168 
170 {}
171 
172 
173 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174 
177 (
178  const label bandI
179 ) const
180 {
181  const fluidMulticomponentThermo& mcThermo =
182  dynamic_cast<const fluidMulticomponentThermo&>(thermo_);
183 
184  const volScalarField& T = thermo_.T();
185  const volScalarField& p = thermo_.p();
186 
187 
189  (
191  (
192  "aCont" + name(bandI),
193  mesh(),
195  extrapolatedCalculatedFvPatchVectorField::typeName
196  )
197  );
198 
199  scalarField& a = ta.ref().primitiveFieldRef();
200 
201  const unitConversion& unitAtm = units()["atm"];
202 
203  forAll(a, celli)
204  {
205  forAllConstIter(HashTable<label>, speciesNames_, iter)
206  {
207  label n = iter();
208  scalar Xipi = 0.0;
209  if (specieIndex_[n] != 0)
210  {
211  // Specie found in the lookUpTable.
212  const volScalarField& ft =
213  mesh_.lookupObject<volScalarField>("ft");
214 
215  const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
216 
217  // moles x pressure [atm]
218  Xipi = unitAtm.toUser(Ynft[specieIndex_[n]]*p[celli]);
219  }
220  else
221  {
222  scalar invWt = 0.0;
223  forAll(mcThermo.Y(), s)
224  {
225  invWt += mcThermo.Y(s)[celli]/mcThermo.WiValue(s);
226  }
227 
228  const label index = mcThermo.species()[iter.key()];
229 
230  const scalar Xk =
231  mcThermo.Y(index)[celli]/(mcThermo.WiValue(index)*invWt);
232 
233  Xipi = unitAtm.toUser(Xk*p[celli]);
234  }
235 
236  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
237 
238  scalar Ti = T[celli];
239  // negative temperature exponents
240  if (coeffs_[n].invTemp())
241  {
242  Ti = 1.0/T[celli];
243  }
244  a[celli] +=
245  Xipi
246  *(
247  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
248  + b[0]
249  );
250  }
251  }
252  ta.ref().correctBoundaryConditions();
253  return ta;
254 }
255 
256 
259 (
260  const label bandI
261 ) const
262 {
263  return aCont(bandI);
264 }
265 
266 
269 (
270  const label bandI
271 ) const
272 {
273  return absorptionEmissionModel::ECont(bandI);
274 }
275 
276 
277 // ************************************************************************* //
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
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 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:139
UPtrList< volScalarField > Yj_
Pointer list of species in the registry involved in the absorption.
Definition: greyMean.H:142
HashTable< label > speciesNames_
Hash table of species names.
Definition: greyMean.H:130
autoPtr< interpolationLookUpTable > lookUpTablePtr_
Look-up table of species related to ft.
Definition: greyMean.H:136
FixedList< label, nSpecies_ > specieIndex_
Indices of species in the look-up table.
Definition: greyMean.H:133
tmp< volScalarField > eCont(const label bandI=0) const
Emission coefficient for continuous phase.
Definition: greyMean.C:259
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMean.C:177
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
Definition: greyMean.C:269
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:197
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
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:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
const dimensionSet dimLength
const HashTable< unitConversion > & units()
Get the table of unit conversions.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p