BilgerMixtureFraction.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) 2020 Thorsten Zirwes
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "basicThermo.H"
31#include "reactingMixture.H"
32#include "thermoPhysicsTypes.H"
33#include "scalarRange.H"
34#include "basicChemistryModel.H"
35#include "psiReactionThermo.H"
36#include "rhoReactionThermo.H"
37#include "BasicChemistryModel.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace functionObjects
45{
48 (
52 );
53}
54}
55
56
57// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58
59void Foam::functionObjects::BilgerMixtureFraction::calcBilgerMixtureFraction()
60{
61 if (!mesh_.foundObject<volScalarField>(resultName_, false))
62 {
64 (
65 IOobject
66 (
67 resultName_,
69 mesh_,
72 ),
73 mesh_,
75 );
76 mesh_.objectRegistry::store(tCo.ptr());
77 }
78
79 auto& f_Bilger = mesh_.lookupObjectRef<volScalarField>(resultName_);
80
81 auto& Y = thermo_.Y();
82
83 f_Bilger = -o2RequiredOx_;
84 forAll(Y, i)
85 {
86 f_Bilger +=
87 Y[i]
88 *(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i] - 0.5*nAtomsO_[i])
89 /thermo_.W(i);
90 }
91 f_Bilger /= o2RequiredFuelOx_;
92 f_Bilger.clip
93 (
96 );
97}
98
99
100bool Foam::functionObjects::BilgerMixtureFraction::readComposition
101(
102 const dictionary& subDict,
103 scalarField& comp
104)
105{
106 auto& Y = thermo_.Y();
107 const speciesTable& speciesTab = thermo_.species();
108
109 // Read mass fractions of all species for the oxidiser or fuel
110 forAll(Y, i)
111 {
112 comp[i] =
113 subDict.getCheckOrDefault<scalar>
114 (
115 speciesTab[i],
116 0,
118 );
119 }
120
121 if (sum(comp) < SMALL)
122 {
124 << "No composition is given" << nl
125 << "Valid species are:" << nl
126 << speciesTab
127 << exit(FatalIOError);
128
129 return false;
130 }
131
132 const word fractionBasisType
133 (
134 subDict.getOrDefault<word>("fractionBasis", "mass")
135 );
136
137 if (fractionBasisType == "mass")
138 {
139 // Normalize fractionBasis to the unity
140 comp /= sum(comp);
141 }
142 else if (fractionBasisType == "mole")
143 {
144 // In case the fractionBasis is given in mole fractions,
145 // convert from mole fractions to normalized mass fractions
146 scalar W(0);
147 forAll(comp, i)
148 {
149 comp[i] *= thermo_.W(i);
150 W += comp[i];
151 }
152 comp /= W;
153 }
154 else
155 {
157 << "The given fractionBasis type is invalid" << nl
158 << "Valid fractionBasis types are" << nl
159 << " \"mass\" (default)" << nl
160 << " \"mole\""
161 << exit(FatalIOError);
162
163 return false;
164 }
165
166 return true;
167}
168
169
170Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Required
171(
172 const scalarField& comp
173) const
174{
175 scalar o2req(0);
176 forAll(thermo_.Y(), i)
177 {
178 o2req +=
179 comp[i]/thermo_.W(i)*(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i]);
180 }
181
182 return o2req;
183}
184
185
186Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Present
187(
188 const scalarField& comp
189) const
190{
191 scalar o2pres(0);
192 forAll(thermo_.Y(), i)
193 {
194 o2pres += comp[i]/thermo_.W(i)*nAtomsO_[i];
195 }
196
197 return 0.5*o2pres;
198}
199
200
201// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202
204(
205 const word& name,
206 const Time& runTime,
207 const dictionary& dict
208)
209:
211 phaseName_(dict.getOrDefault<word>("phase", word::null)),
212 resultName_
213 (
214 dict.getOrDefault<word>
215 (
216 "result",
217 IOobject::groupName("f_Bilger", phaseName_)
218 )
219 ),
220 thermo_
221 (
222 mesh_.lookupObject<basicSpecieMixture>
223 (
224 IOobject::groupName(basicThermo::dictName, phaseName_)
225 )
226 ),
227 nSpecies_(thermo_.Y().size()),
228 o2RequiredOx_(0),
229 o2RequiredFuelOx_(0),
230 nAtomsC_(nSpecies_, 0),
231 nAtomsS_(nSpecies_, 0),
232 nAtomsH_(nSpecies_, 0),
233 nAtomsO_(nSpecies_, 0),
234 Yoxidiser_(nSpecies_, 0),
235 Yfuel_(nSpecies_, 0)
236{
237 read(dict);
238
239 calcBilgerMixtureFraction();
240}
241
242
243// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244
246{
248 {
249 return false;
250 }
251
252 Info<< nl << type() << " " << name() << ":" << nl;
253
254 phaseName_ = dict.getOrDefault<word>("phase", word::null);
255 resultName_ =
256 dict.getOrDefault<word>
257 (
258 "result",
259 IOobject::groupName("f_Bilger", phaseName_)
260 );
261
262 nSpecies_ = thermo_.Y().size();
263
264 if (nSpecies_ == 0)
265 {
267 << "Number of input species is zero"
268 << exit(FatalError);
269 }
270
271 nAtomsC_.resize(nSpecies_, 0);
272 nAtomsS_.resize(nSpecies_, 0);
273 nAtomsH_.resize(nSpecies_, 0);
274 nAtomsO_.resize(nSpecies_, 0);
275
276 auto& Y = thermo_.Y();
277 const speciesTable& speciesTab = thermo_.species();
278
279 typedef BasicChemistryModel<psiReactionThermo> psiChemistryModelType;
280 typedef BasicChemistryModel<rhoReactionThermo> rhoChemistryModelType;
281
282 const auto* psiChemPtr =
283 mesh_.cfindObject<psiChemistryModelType>("chemistryProperties");
284
285 const auto* rhoChemPtr =
286 mesh_.cfindObject<rhoChemistryModelType>("chemistryProperties");
287
289
290 if (psiChemPtr)
291 {
292 speciesCompPtr.reset((*psiChemPtr).thermo().specieComposition());
293 }
294 else if (rhoChemPtr)
295 {
296 speciesCompPtr.reset((*rhoChemPtr).thermo().specieComposition());
297 }
298 else
299 {
301 << "BasicChemistryModel not found"
302 << exit(FatalError);
303 }
304
305 forAll(Y, i)
306 {
307 const List<specieElement>& curSpecieComposition =
308 (speciesCompPtr.ref())[speciesTab[i]];
309
310 forAll(curSpecieComposition, j)
311 {
312 const word& e = curSpecieComposition[j].name();
313 const label nAtoms = curSpecieComposition[j].nAtoms();
314
315 if (e == "C")
316 {
317 nAtomsC_[i] = nAtoms;
318 }
319 else if (e == "S")
320 {
321 nAtomsS_[i] = nAtoms;
322 }
323 else if (e == "H")
324 {
325 nAtomsH_[i] = nAtoms;
326 }
327 else if (e == "O")
328 {
329 nAtomsO_[i] = nAtoms;
330 }
331 }
332 }
333
334 if (sum(nAtomsO_) == 0)
335 {
337 << "No specie contains oxygen"
338 << exit(FatalError);
339 }
340
341 Yoxidiser_.resize(nSpecies_, 0);
342 Yfuel_.resize(nSpecies_, 0);
343
344 if
345 (
346 !readComposition(dict.subDict("oxidiser"), Yoxidiser_)
347 || !readComposition(dict.subDict("fuel"), Yfuel_)
348 )
349 {
350 return false;
351 }
352
353 o2RequiredOx_ = o2Required(Yoxidiser_) - o2Present(Yoxidiser_);
354
355 if (o2RequiredOx_ > 0)
356 {
358 << "Oxidiser composition contains not enough oxygen" << endl
359 << "Mixed up fuel and oxidiser compositions?"
360 << exit(FatalError);
361 }
362
363 const scalar o2RequiredFuel = o2Required(Yfuel_) - o2Present(Yfuel_);
364
365 if (o2RequiredFuel < 0)
366 {
368 << "Fuel composition contains too much oxygen" << endl
369 << "Mixed up fuel and oxidiser compositions?"
370 << exit(FatalError);
371 }
372
373 o2RequiredFuelOx_ = o2RequiredFuel - o2RequiredOx_;
374
375 if (mag(o2RequiredFuelOx_) < SMALL)
376 {
378 << "Fuel and oxidiser have the same composition"
379 << exit(FatalError);
380 }
381
382 return true;
383}
384
385
387{
388 calcBilgerMixtureFraction();
389
390 return true;
391}
392
393
395{
396 return clearObject(resultName_);
397}
398
399
401{
402 Log << type() << " " << name() << " write:" << nl
403 << " writing field " << resultName_ << endl;
404
405 return writeObject(resultName_);
406}
407
408
409// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Basic chemistry model templated on thermodynamics.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Calculates the Bilger mixture-fraction field (i.e. ) based on the elemental composition of the mixtur...
virtual bool clear()
Clear the Bilger mixture-fraction field from registry.
virtual bool execute()
Calculate the Bilger mixture-fraction field.
virtual bool write()
Write Bilger mixture-fraction field.
virtual bool read(const dictionary &)
Read the BilgerMixtureFraction data.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
Definition: mag.H:153
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Type & lookupObjectRef(const word &name, const bool recursive=false) const
static constexpr scalarRange ge0() noexcept
A greater-equals zero bound.
Definition: scalarRangeI.H:89
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
PtrList< volScalarField > & Y
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const word dictName("faMeshDefinition")
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
hashedWordList speciesTable
A table of species as a hashedWordList.
Definition: speciesTable.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Type definitions for thermo-physics models.