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