veryInhomogeneousMixture.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 {
34  "ft",
35  "fu",
36  "b"
37 };
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class ThermoType>
44 (
45  const dictionary& thermoDict,
46  const fvMesh& mesh,
47  const word& phaseName
48 )
49 :
51  (
52  thermoDict,
53  speciesTable(nSpecies_, specieNames_),
54  mesh,
55  phaseName
56  ),
57 
58  stoicRatio_(thermoDict.lookup("stoichiometricAirFuelMassRatio")),
59 
60  fuel_(thermoDict.subDict("fuel")),
61  oxidant_(thermoDict.subDict("oxidant")),
62  products_(thermoDict.subDict("burntProducts")),
63 
64  mixture_("mixture", fuel_),
65 
66  ft_(Y("ft")),
67  fu_(Y("fu")),
68  b_(Y("b"))
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class ThermoType>
76 (
77  const scalar ft,
78  const scalar fu
79 ) const
80 {
81  if (ft < 0.0001)
82  {
83  return oxidant_;
84  }
85  else
86  {
87  scalar ox = 1 - ft - (ft - fu)*stoicRatio().value();
88  scalar pr = 1 - fu - ox;
89 
90  mixture_ = fu/fuel_.W()*fuel_;
91  mixture_ += ox/oxidant_.W()*oxidant_;
92  mixture_ += pr/products_.W()*products_;
93 
94  return mixture_;
95  }
96 }
97 
98 
99 template<class ThermoType>
101 (
102  const dictionary& thermoDict
103 )
104 {
105  fuel_ = ThermoType(thermoDict.subDict("fuel"));
106  oxidant_ = ThermoType(thermoDict.subDict("oxidant"));
107  products_ = ThermoType(thermoDict.subDict("burntProducts"));
108 }
109 
110 
111 template<class ThermoType>
113 (
114  const label speciei
115 ) const
116 {
117  if (speciei == 0)
118  {
119  return fuel_;
120  }
121  else if (speciei == 1)
122  {
123  return oxidant_;
124  }
125  else if (speciei == 2)
126  {
127  return products_;
128  }
129  else
130  {
132  << "Unknown specie index " << speciei << ". Valid indices are 0..2"
133  << abort(FatalError);
134 
135  return fuel_;
136  }
137 }
138 
139 
140 // ************************************************************************* //
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
hashedWordList speciesTable
A table of species as a hashedWordList.
Definition: speciesTable.H:41
void read(const dictionary &)
Read dictionary.
A class for handling words, derived from string.
Definition: word.H:59
Specialization of the basicSpecieMixture for combustion.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
const ThermoType & mixture(const scalar, const scalar) const
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::veryInhomogeneousMixture.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451