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-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 "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 scalar ft,
50  const scalar b
51 ) const
52 {
53  if (ft < 0.0001)
54  {
55  return oxidant_;
56  }
57  else
58  {
59  const scalar fu = b*ft + (1 - b)*fres(ft);
60  const scalar ox = 1 - ft - (ft - fu)*stoicRatio_.value();
61  const scalar pr = 1 - fu - ox;
62 
63  mixture_ = fu*fuel_;
64  mixture_ += ox*oxidant_;
65  mixture_ += pr*products_;
66 
67  return mixture_;
68  }
69 }
70 
71 
72 template<class ThermoType>
75 (
76  const scalarFieldListSlice& Y
77 ) const
78 {
79  return mixture(Y[FT], Y[B]);
80 }
81 
82 
83 template<class ThermoType>
86 (
87  const scalarFieldListSlice& Y
88 ) const
89 {
90  return mixture(Y[FT], Y[B]);
91 }
92 
93 
94 template<class ThermoType>
97 (
98  const scalarFieldListSlice&,
100 ) const
101 {
102  return mixture;
103 }
104 
105 
106 template<class ThermoType>
109 (
110  const scalarFieldListSlice& Y
111 ) const
112 {
113  return mixture(Y[FT], 1);
114 }
115 
116 
117 template<class ThermoType>
120 (
121  const scalarFieldListSlice& Y
122 ) const
123 {
124  return mixture(Y[FT], 0);
125 }
126 
127 
128 template<class ThermoType>
130 {
131  stoicRatio_ =
132  dimensionedScalar("stoichiometricAirFuelMassRatio", dimless, dict);
133 
134  fuel_ = ThermoType("fuel", dict.subDict("fuel"));
135  oxidant_ = ThermoType("oxidant", dict.subDict("oxidant"));
136  products_ = ThermoType("burntProducts", dict.subDict("burntProducts"));
137 }
138 
139 
140 // ************************************************************************* //
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 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.
const thermoType & products(const scalarFieldListSlice &) const
Return the product mixture.
ThermoType transportMixtureType
Mixing type for transport properties.
inhomogeneousMixture(const dictionary &)
Construct from a dictionary.
ThermoType thermoMixtureType
Mixing type for thermodynamic properties.
const thermoType & reactants(const scalarFieldListSlice &) const
Return the reactant mixture.
const transportMixtureType & transportMixture(const scalarFieldListSlice &) const
Return the mixture for transport properties.
void read(const dictionary &)
Read dictionary.
const thermoMixtureType & thermoMixture(const scalarFieldListSlice &) const
Return the mixture for thermodynamic properties.
volScalarField & b
Definition: createFields.H:25
const dimensionSet dimless
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
PtrList< volScalarField > & Y