Loading...
Searching...
No Matches
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-2024 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"
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 auto* resultPtr = mesh_.getObjectPtr<volScalarField>(resultName_);
62
63 if (!resultPtr)
64 {
65 resultPtr = new volScalarField
66 (
68 (
69 resultName_,
71 mesh_.thisDb(),
75 ),
76 mesh_,
78 );
79 regIOobject::store(resultPtr);
80 }
81 auto& f_Bilger = *resultPtr;
82
83 auto& Y = thermo_.Y();
84
85 f_Bilger = -o2RequiredOx_;
86 forAll(Y, i)
87 {
88 f_Bilger +=
89 Y[i]
90 *(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i] - 0.5*nAtomsO_[i])
91 /thermo_.W(i);
92 }
93 f_Bilger /= o2RequiredFuelOx_;
94 f_Bilger.clamp_range(zero_one{});
95}
96
97
98bool Foam::functionObjects::BilgerMixtureFraction::readComposition
99(
100 const dictionary& subDict,
101 scalarField& comp
102)
103{
104 auto& Y = thermo_.Y();
105 const speciesTable& speciesTab = thermo_.species();
106
107 // Read mass fractions of all species for the oxidiser or fuel
108 forAll(Y, i)
109 {
110 comp[i] =
111 subDict.getCheckOrDefault<scalar>
112 (
113 speciesTab[i],
114 0,
116 );
117 }
118
119 if (sum(comp) < SMALL)
120 {
122 << "No composition is given" << nl
123 << "Valid species are:" << nl
124 << speciesTab
125 << exit(FatalIOError);
126
127 return false;
128 }
129
130 const word fractionBasisType
131 (
132 subDict.getOrDefault<word>("fractionBasis", "mass")
133 );
134
135 if (fractionBasisType == "mass")
136 {
137 // Normalize fractionBasis to the unity
138 comp /= sum(comp);
139 }
140 else if (fractionBasisType == "mole")
141 {
142 // In case the fractionBasis is given in mole fractions,
143 // convert from mole fractions to normalized mass fractions
144 scalar W(0);
145 forAll(comp, i)
146 {
147 comp[i] *= thermo_.W(i);
148 W += comp[i];
149 }
150 comp /= W;
151 }
152 else
153 {
155 << "The given fractionBasis type is invalid" << nl
156 << "Valid fractionBasis types are" << nl
157 << " \"mass\" (default)" << nl
158 << " \"mole\""
159 << exit(FatalIOError);
160
161 return false;
162 }
163
164 return true;
165}
166
167
168Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Required
169(
170 const scalarField& comp
171) const
172{
173 scalar o2req(0);
174 forAll(thermo_.Y(), i)
175 {
176 o2req +=
177 comp[i]/thermo_.W(i)*(nAtomsC_[i] + nAtomsS_[i] + 0.25*nAtomsH_[i]);
178 }
179
180 return o2req;
181}
182
183
184Foam::scalar Foam::functionObjects::BilgerMixtureFraction::o2Present
185(
186 const scalarField& comp
187) const
188{
189 scalar o2pres(0);
190 forAll(thermo_.Y(), i)
191 {
192 o2pres += comp[i]/thermo_.W(i)*nAtomsO_[i];
193 }
195 return 0.5*o2pres;
196}
197
198
199// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
200
202(
203 const word& name,
204 const Time& runTime,
205 const dictionary& dict
206)
207:
209 phaseName_(dict.getOrDefault<word>("phase", word::null)),
210 resultName_
211 (
212 dict.getOrDefault<word>
213 (
214 "result",
215 IOobject::groupName("f_Bilger", phaseName_)
216 )
217 ),
218 thermo_
219 (
220 mesh_.lookupObject<basicSpecieMixture>
221 (
222 IOobject::groupName(basicThermo::dictName, phaseName_)
223 )
224 ),
225 nSpecies_(thermo_.Y().size()),
226 o2RequiredOx_(0),
227 o2RequiredFuelOx_(0),
228 nAtomsC_(nSpecies_, 0),
229 nAtomsS_(nSpecies_, 0),
230 nAtomsH_(nSpecies_, 0),
231 nAtomsO_(nSpecies_, 0),
232 Yoxidiser_(nSpecies_, 0),
233 Yfuel_(nSpecies_, 0)
234{
235 read(dict);
237 calcBilgerMixtureFraction();
238}
239
240
241// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242
244{
246 {
247 return false;
248 }
249
250 Info<< nl << type() << " " << name() << ":" << nl;
251
252 phaseName_ = dict.getOrDefault<word>("phase", word::null);
253 resultName_ =
254 dict.getOrDefault<word>
255 (
256 "result",
257 IOobject::groupName("f_Bilger", phaseName_)
258 );
259
260 nSpecies_ = thermo_.Y().size();
261
262 if (nSpecies_ == 0)
263 {
265 << "Number of input species is zero"
266 << exit(FatalError);
267 }
268
269 nAtomsC_.resize(nSpecies_, 0);
270 nAtomsS_.resize(nSpecies_, 0);
271 nAtomsH_.resize(nSpecies_, 0);
272 nAtomsO_.resize(nSpecies_, 0);
273
274 auto& Y = thermo_.Y();
275 const speciesTable& speciesTab = thermo_.species();
276
277 typedef BasicChemistryModel<psiReactionThermo> psiChemistryModelType;
278 typedef BasicChemistryModel<rhoReactionThermo> rhoChemistryModelType;
279
280 const auto* psiChemPtr =
281 mesh_.cfindObject<psiChemistryModelType>("chemistryProperties");
282
283 const auto* rhoChemPtr =
284 mesh_.cfindObject<rhoChemistryModelType>("chemistryProperties");
285
286 autoPtr<speciesCompositionTable> speciesCompPtr;
287
288 if (psiChemPtr)
289 {
290 speciesCompPtr.reset((*psiChemPtr).thermo().specieComposition());
291 }
292 else if (rhoChemPtr)
293 {
294 speciesCompPtr.reset((*rhoChemPtr).thermo().specieComposition());
295 }
296 else
297 {
299 << "BasicChemistryModel not found"
300 << exit(FatalError);
301 }
302
303 forAll(Y, i)
304 {
305 const List<specieElement>& curSpecieComposition =
306 (speciesCompPtr.ref())[speciesTab[i]];
307
308 forAll(curSpecieComposition, j)
309 {
310 const word& e = curSpecieComposition[j].name();
311 const label nAtoms = curSpecieComposition[j].nAtoms();
312
313 if (e == "C")
314 {
315 nAtomsC_[i] = nAtoms;
316 }
317 else if (e == "S")
318 {
319 nAtomsS_[i] = nAtoms;
320 }
321 else if (e == "H")
322 {
323 nAtomsH_[i] = nAtoms;
324 }
325 else if (e == "O")
326 {
327 nAtomsO_[i] = nAtoms;
328 }
329 }
330 }
331
332 if (sum(nAtomsO_) == 0)
333 {
335 << "No specie contains oxygen"
336 << exit(FatalError);
337 }
338
339 Yoxidiser_.resize(nSpecies_, 0);
340 Yfuel_.resize(nSpecies_, 0);
341
342 if
343 (
344 !readComposition(dict.subDict("oxidiser"), Yoxidiser_)
345 || !readComposition(dict.subDict("fuel"), Yfuel_)
346 )
347 {
348 return false;
349 }
350
351 o2RequiredOx_ = o2Required(Yoxidiser_) - o2Present(Yoxidiser_);
352
353 if (o2RequiredOx_ > 0)
354 {
356 << "Oxidiser composition contains not enough oxygen" << endl
357 << "Mixed up fuel and oxidiser compositions?"
358 << exit(FatalError);
359 }
360
361 const scalar o2RequiredFuel = o2Required(Yfuel_) - o2Present(Yfuel_);
362
363 if (o2RequiredFuel < 0)
364 {
366 << "Fuel composition contains too much oxygen" << endl
367 << "Mixed up fuel and oxidiser compositions?"
368 << exit(FatalError);
369 }
370
371 o2RequiredFuelOx_ = o2RequiredFuel - o2RequiredOx_;
372
373 if (mag(o2RequiredFuelOx_) < SMALL)
374 {
376 << "Fuel and oxidiser have the same composition"
378 }
379
380 return true;
381}
382
383
386 calcBilgerMixtureFraction();
387
388 return true;
389}
390
393{
394 return clearObject(resultName_);
395}
396
397
399{
400 Log << type() << " " << name() << " write:" << nl
401 << " writing field " << resultName_ << endl;
402
403 return writeObject(resultName_);
404}
405
406
407// ************************************************************************* //
#define Log
Definition PDRblock.C:28
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.
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Calculates the Bilger mixture-fraction field (i.e. ) based on the elemental composition of the mixtur...
BilgerMixtureFraction(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
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.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool writeObject(const word &fieldName)
Write field if present in the (sub) objectRegistry.
bool clearObject(const word &fieldName)
Clear field from the (sub) objectRegistry if present.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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...
bool store()
Register object with its registry and transfer ownership to the registry.
static constexpr scalarRange ge0() noexcept
A greater-equals zero bound.
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
PtrList< volScalarField > & Y
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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 POSIX.C:801
hashedWordList speciesTable
A table of species as a hashedWordList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Type definitions for thermo-physics models.