leanInhomogeneousMixture.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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ThermoType>
32 (
33  const dictionary& dict
34 )
35 :
36  stoicRatio_("stoichiometricAirFuelMassRatio", dimless, dict),
37  fuel_("fuel", dict.subDict("fuel")),
38  oxidant_("oxidant", dict.subDict("oxidant")),
39  products_("burntProducts", dict.subDict("burntProducts")),
40  mixture_("mixture", fuel_)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
46 template<class ThermoType>
48 (
49  const scalar ft
50 ) const
51 {
52  return max(ft - (scalar(1) - ft)/stoicRatio_.value(), scalar(0));
53 }
54 
55 
56 template<class ThermoType>
58 (
59  const scalarFieldListSlice& Y
60 ) const
61 {
62  return fres(Y[FT]);
63 }
64 
65 
66 template<class ThermoType>
68 (
69  const scalar ft,
70  const scalar b
71 ) const
72 {
73  if (ft < 0.0001)
74  {
75  return oxidant_;
76  }
77  else
78  {
79  const scalar fu = b*ft + (1 - b)*fres(ft);
80  const scalar ox = 1 - ft - (ft - fu)*stoicRatio_.value();
81  const scalar pr = 1 - fu - ox;
82 
83  mixture_ = fu*fuel_;
84  mixture_ += ox*oxidant_;
85  mixture_ += pr*products_;
86 
87  return mixture_;
88  }
89 }
90 
91 
92 template<class ThermoType>
95 (
96  const scalarFieldListSlice& Y
97 ) const
98 {
99  return mixture(Y[FT], Y[B]);
100 }
101 
102 
103 template<class ThermoType>
106 (
107  const scalarFieldListSlice& Y
108 ) const
109 {
110  return mixture(Y[FT], Y[B]);
111 }
112 
113 
114 template<class ThermoType>
117 (
118  const scalarFieldListSlice&,
120 ) const
121 {
122  return mixture;
123 }
124 
125 
126 template<class ThermoType>
129 (
130  const scalarFieldListSlice& Y
131 ) const
132 {
133  return mixture(Y[FT], 1);
134 }
135 
136 
137 template<class ThermoType>
140 (
141  const scalarFieldListSlice& Y
142 ) const
143 {
144  return mixture(Y[FT], 0);
145 }
146 
147 
148 template<class ThermoType>
150 {
151  stoicRatio_ =
152  dimensionedScalar("stoichiometricAirFuelMassRatio", dimless, dict);
153 
154  fuel_ = ThermoType("fuel", dict.subDict("fuel"));
155  oxidant_ = ThermoType("oxidant", dict.subDict("oxidant"));
156  products_ = ThermoType("burntProducts", dict.subDict("burntProducts"));
157 }
158 
159 
160 // ************************************************************************* //
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
const thermoType & mixture(const scalar ft, const scalar b) const
Return the mixture for the given composition.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
leanInhomogeneousMixture(const dictionary &)
Construct from a dictionary.
ThermoType transportMixtureType
Mixing type for transport properties.
ThermoType thermoMixtureType
Mixing type for thermodynamic properties.
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
const thermoType & products(const scalarFieldListSlice &) const
Return the product mixture.
void read(const dictionary &)
Read dictionary.
scalar fres(const scalar ft) const
Return the residual fraction of fuel in the burnt mixture.
const thermoType & reactants(const scalarFieldListSlice &) const
Return the reactant mixture.
volScalarField & b
Definition: createFields.H:27
static const coefficient B("B", dimless, 18.678)
const dimensionSet dimless
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
PtrList< volScalarField > & Y