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-2020 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 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 {
34  "ft",
35  "b"
36 };
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class ThermoType>
43 (
44  const dictionary& thermoDict,
45  const fvMesh& mesh,
46  const word& phaseName
47 )
48 :
50  (
51  thermoDict,
52  speciesTable(nSpecies_, specieNames_),
53  mesh,
54  phaseName
55  ),
56 
57  stoicRatio_(thermoDict.lookup("stoichiometricAirFuelMassRatio")),
58 
59  fuel_(thermoDict.subDict("fuel")),
60  oxidant_(thermoDict.subDict("oxidant")),
61  products_(thermoDict.subDict("burntProducts")),
62 
63  mixture_("mixture", fuel_),
64 
65  ft_(Y("ft")),
66  b_(Y("b"))
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 template<class ThermoType>
74 (
75  const scalar ft,
76  const scalar b
77 ) const
78 {
79  if (ft < 0.0001)
80  {
81  return oxidant_;
82  }
83  else
84  {
85  scalar fu = b*ft + (1.0 - b)*fres(ft, stoicRatio().value());
86  scalar ox = 1 - ft - (ft - fu)*stoicRatio().value();
87  scalar pr = 1 - fu - ox;
88 
89  mixture_ = fu*fuel_;
90  mixture_ += ox*oxidant_;
91  mixture_ += pr*products_;
92 
93  return mixture_;
94  }
95 }
96 
97 
98 template<class ThermoType>
100 {
101  stoicRatio_ = thermoDict.lookup("stoichiometricAirFuelMassRatio");
102 
103  fuel_ = ThermoType(thermoDict.subDict("fuel"));
104  oxidant_ = ThermoType(thermoDict.subDict("oxidant"));
105  products_ = ThermoType(thermoDict.subDict("burntProducts"));
106 }
107 
108 
109 template<class ThermoType>
111 (
112  const label speciei
113 ) const
114 {
115  if (speciei == 0)
116  {
117  return fuel_;
118  }
119  else if (speciei == 1)
120  {
121  return oxidant_;
122  }
123  else if (speciei == 2)
124  {
125  return products_;
126  }
127  else
128  {
130  << "Unknown specie index " << speciei << ". Valid indices are 0..2"
131  << abort(FatalError);
132 
133  return fuel_;
134  }
135 }
136 
137 
138 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const ThermoType & specieThermo(const label speciei) const
Return thermo based on index.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Foam::inhomogeneousMixture.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
hashedWordList speciesTable
A table of species as a hashedWordList.
Definition: speciesTable.H:41
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
inhomogeneousMixture(const dictionary &, const fvMesh &, const word &)
Construct from dictionary, mesh and phase name.
A class for handling words, derived from string.
Definition: word.H:59
Specialisation of the basicMixture for combustion.
void read(const dictionary &)
Read dictionary.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const thermoType & mixture(const scalar, const scalar) const
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844