egrMixture.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-2023 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 "egrMixture.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ThermoType>
32 :
33  stoicRatio_("stoichiometricAirFuelMassRatio", dimless, dict),
34  fuel_("fuel", dict.subDict("fuel")),
35  oxidant_("oxidant", dict.subDict("oxidant")),
36  products_("burntProducts", dict.subDict("burntProducts")),
37  mixture_("mixture", fuel_)
38 {}
39 
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
43 template<class ThermoType>
45 (
46  const scalar ft,
47  const scalar b,
48  const scalar egr
49 ) const
50 {
51  if (ft < 0.0001 && egr < 0.0001)
52  {
53  return oxidant_;
54  }
55  else
56  {
57  scalar fu = b*ft + (1 - b)*fres(ft);
58  scalar ox = 1 - ft - (ft - fu)*stoicRatio_.value();
59 
60  fu *= 1 - egr;
61  ox *= 1 - egr;
62 
63  const scalar pr = 1 - fu - ox;
64 
65  mixture_ = fu*fuel_;
66  mixture_ += ox*oxidant_;
67  mixture_ += pr*products_;
68 
69  return mixture_;
70  }
71 }
72 
73 
74 template<class ThermoType>
77 (
78  const scalarFieldListSlice& Y
79 ) const
80 {
81  return mixture(Y[FT], Y[B], Y[EGR]);
82 }
83 
84 
85 template<class ThermoType>
88 (
89  const scalarFieldListSlice& Y
90 ) const
91 {
92  return mixture(Y[FT], Y[B], Y[EGR]);
93 }
94 
95 
96 template<class ThermoType>
99 (
100  const scalarFieldListSlice&,
102 ) const
103 {
104  return mixture;
105 }
106 
107 
108 template<class ThermoType>
110 Foam::egrMixture<ThermoType>::reactants(const scalarFieldListSlice& Y) const
111 {
112  return mixture(Y[FT], 1, Y[EGR]);
113 }
114 
115 
116 template<class ThermoType>
118 Foam::egrMixture<ThermoType>::products(const scalarFieldListSlice& Y) const
119 {
120  return mixture(Y[FT], 0, 0);
121 }
122 
123 
124 template<class ThermoType>
126 {
127  stoicRatio_ =
128  dimensionedScalar("stoichiometricAirFuelMassRatio", dimless, dict);
129 
130  fuel_ = ThermoType("fuel", dict.subDict("fuel"));
131  oxidant_ = ThermoType("oxidant", dict.subDict("oxidant"));
132  products_ = ThermoType("burntProducts", dict.subDict("burntProducts"));
133 }
134 
135 
136 // ************************************************************************* //
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
Definition: egrMixture.C:88
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: egrMixture.H:57
egrMixture(const dictionary &)
Construct from a dictionary.
Definition: egrMixture.C:31
ThermoType transportMixtureType
Mixing type for transport properties.
Definition: egrMixture.H:63
ThermoType thermoMixtureType
Mixing type for thermodynamic properties.
Definition: egrMixture.H:60
const thermoType & mixture(const scalar ft, const scalar b, const scalar egr) const
Return the mixture for the given composition.
Definition: egrMixture.C:45
void read(const dictionary &)
Read dictionary.
Definition: egrMixture.C:125
const thermoType & products(const scalarFieldListSlice &) const
Return the product mixture.
Definition: egrMixture.C:118
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
Definition: egrMixture.C:77
const thermoType & reactants(const scalarFieldListSlice &) const
Return the reactant mixture.
Definition: egrMixture.C:110
volScalarField & b
Definition: createFields.H:25
const dimensionSet dimless
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
PtrList< volScalarField > & Y