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-2021 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"
28 #include "unitConversion.H"
30 #include "basicSpecieMixture.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace radiationModels
37 {
38 namespace absorptionEmissionModels
39 {
41 
43  (
45  greyMean,
47  );
48 }
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const dictionary& dict,
58  const fvMesh& mesh,
59  const word& modelName
60 )
61 :
63  coeffsDict_(dict.subDict(modelName + "Coeffs")),
64  speciesNames_(0),
65  specieIndex_(label(0)),
66  lookUpTablePtr_(),
67  thermo_(mesh.lookupObject<fluidThermo>(physicalProperties::typeName)),
68  Yj_(nSpecies_)
69 {
70  if (!isA<basicSpecieMixture>(thermo_))
71  {
73  << "Model requires a multi-component thermo package"
74  << abort(FatalError);
75  }
76 
77  label nFunc = 0;
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;
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 label bandI
181 ) const
182 {
183  const basicSpecieMixture& mixture =
184  dynamic_cast<const basicSpecieMixture&>(thermo_);
185 
186  const volScalarField& T = thermo_.T();
187  const volScalarField& p = thermo_.p();
188 
189 
191  (
193  (
194  "aCont" + name(bandI),
195  mesh(),
197  extrapolatedCalculatedFvPatchVectorField::typeName
198  )
199  );
200 
201  scalarField& a = ta.ref().primitiveFieldRef();
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  // moles x pressure [atm]
217  Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
218  }
219  else
220  {
221  scalar invWt = 0.0;
222  forAll(mixture.Y(), s)
223  {
224  invWt += mixture.Y(s)[celli]/mixture.Wi(s);
225  }
226 
227  label index = mixture.species()[iter.key()];
228  scalar Xk = mixture.Y(index)[celli]/(mixture.Wi(index)*invWt);
229 
230  Xipi = Xk*paToAtm(p[celli]);
231  }
232 
233  const absorptionCoeffs::coeffArray& b = coeffs_[n].coeffs(T[celli]);
234 
235  scalar Ti = T[celli];
236  // negative temperature exponents
237  if (coeffs_[n].invTemp())
238  {
239  Ti = 1.0/T[celli];
240  }
241  a[celli] +=
242  Xipi
243  *(
244  ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
245  + b[0]
246  );
247  }
248  }
249  ta.ref().correctBoundaryConditions();
250  return ta;
251 }
252 
253 
256 (
257  const label bandI
258 ) const
259 {
260  return aCont(bandI);
261 }
262 
263 
266 (
267  const label bandI
268 ) const
269 {
270  return absorptionEmissionModel::ECont(bandI);
271 }
272 
273 
274 // ************************************************************************* //
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 >> &)
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
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
A class for handling file names.
Definition: fileName.H:82
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
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.
An abstract 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:256
tmp< volScalarField > aCont(const label bandI=0) const
Absorption coefficient for continuous phase.
Definition: greyMean.C:179
tmp< volScalarField > ECont(const label bandI=0) const
Emission contribution for continuous phase.
Definition: greyMean.C:266
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:56
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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:27
addToRunTimeSelectionTable(absorptionEmissionModel, greyMeanCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
error FatalError
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p
Unit conversion functions.