inhomogeneousMixture.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 
26 #include "inhomogeneousMixture.H"
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 scalarFieldListSlice& Y
50 ) const
51 {
52  return max(Y[FT] - (scalar(1) - Y[FT])/stoicRatio_.value(), scalar(0));
53 }
54 
55 
56 template<class ThermoType>
58 (
59  const scalar ft,
60  const scalar fu
61 ) const
62 {
63  if (ft < 0.0001)
64  {
65  return oxidant_;
66  }
67  else
68  {
69  const scalar ox = 1 - ft - (ft - fu)*stoicRatio_.value();
70  const scalar pr = 1 - fu - ox;
71 
72  mixture_ = fu*fuel_;
73  mixture_ += ox*oxidant_;
74  mixture_ += pr*products_;
75 
76  return mixture_;
77  }
78 }
79 
80 
81 template<class ThermoType>
84 (
85  const scalarFieldListSlice& Y
86 ) const
87 {
88  return mixture(Y[FT], Y[FU]);
89 }
90 
91 
92 template<class ThermoType>
95 (
96  const scalarFieldListSlice& Y
97 ) const
98 {
99  return mixture(Y[FT], Y[FU]);
100 }
101 
102 
103 template<class ThermoType>
106 (
107  const scalarFieldListSlice&,
109 ) const
110 {
111  return mixture;
112 }
113 
114 
115 template<class ThermoType>
118 (
119  const scalarFieldListSlice& Y
120 ) const
121 {
122  return mixture(Y[FT], Y[FT]);
123 }
124 
125 
126 template<class ThermoType>
129 (
130  const scalarFieldListSlice& Y
131 ) const
132 {
133  return mixture(Y[FT], fres(Y));
134 }
135 
136 
137 template<class ThermoType>
139 (
140  const dictionary& dict
141 )
142 {
143  stoicRatio_ =
144  dimensionedScalar("stoichiometricAirFuelMassRatio", dimless, dict);
145 
146  fuel_ = ThermoType("fuel", dict.subDict("fuel"));
147  oxidant_ = ThermoType("oxidant", dict.subDict("oxidant"));
148  products_ = ThermoType("burntProducts", dict.subDict("burntProducts"));
149 }
150 
151 
152 // ************************************************************************* //
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 &Y) const
Return the mixture for thermodynamic properties.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
const thermoType & products(const scalarFieldListSlice &) const
Return the product mixture.
const transportMixtureType & transportMixture(const scalarFieldListSlice &Y) const
Return the mixture for transport properties.
ThermoType transportMixtureType
Mixing type for transport properties.
scalar fres(const scalarFieldListSlice &Y) const
Return the residual fraction of fuel in the burnt mixture.
inhomogeneousMixture(const dictionary &)
Construct from a dictionary.
ThermoType thermoMixtureType
Mixing type for thermodynamic properties.
const thermoType & reactants(const scalarFieldListSlice &) const
Return the reactant mixture.
void read(const dictionary &)
Read dictionary.
const thermoType & mixture(const scalar ft, const scalar fu) const
Return the mixture for the given composition.
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