Loading...
Searching...
No Matches
egrMixture.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
28#include "egrMixture.H"
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class ThermoType>
34Foam::egrMixture<ThermoType>::egrMixture
35(
37 const fvMesh& mesh,
38 const word& phaseName
39)
40:
42 (
44 speciesTable({"ft", "b", "egr"}),
45 mesh,
46 phaseName
47 ),
48
49 stoicRatio_("stoichiometricAirFuelMassRatio", dimless, thermoDict),
50
51 fuel_(thermoDict.subDict("fuel")),
52 oxidant_(thermoDict.subDict("oxidant")),
53 products_(thermoDict.subDict("burntProducts")),
54
55 mixture_("mixture", fuel_),
56
57 ft_(Y("ft")),
58 b_(Y("b")),
59 egr_(Y("egr"))
60{}
61
62
63// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64
65template<class ThermoType>
67(
68 const scalar ft,
69 const scalar b,
70 const scalar egr
71) const
72{
73 if (ft < 0.0001)
74 {
75 return oxidant_;
76 }
77 else
78 {
79
80 scalar fu = b*ft + (1.0 - b)*fres(ft, stoicRatio().value());
81 scalar ox = 1 - ft - (ft - fu)*stoicRatio().value();
82
83 fu *= (1.0 - egr);
84 ox *= (1.0 - egr);
85
86 scalar pr = 1 - fu - ox;
87
88 mixture_ = fu*fuel_;
89 mixture_ += ox*oxidant_;
90 mixture_ += pr*products_;
92 return mixture_;
93 }
94}
95
96
97template<class ThermoType>
99{
100 thermoDict.readEntry("stoichiometricAirFuelMassRatio", stoicRatio_);
101
102 fuel_ = ThermoType(thermoDict.subDict("fuel"));
103 oxidant_ = ThermoType(thermoDict.subDict("oxidant"));
104 products_ = ThermoType(thermoDict.subDict("burntProducts"));
105}
106
107
108template<class ThermoType>
110(
111 const label speciei
112) const
113{
114 if (speciei == 0)
115 {
116 return fuel_;
117 }
118 else if (speciei == 1)
119 {
120 return oxidant_;
121 }
122 else if (speciei == 2)
123 {
124 return products_;
125 }
126 else
127 {
129 << "Unknown specie index " << speciei << ". Valid indices are 0..2"
130 << abort(FatalError);
131
132 return fuel_;
133 }
134}
135
136
137// ************************************************************************* //
basicCombustionMixture(const dictionary &thermoDict, const wordList &specieNames, const fvMesh &mesh, const word &phaseName)
Construct from dictionary, specie names, mesh and phase name.
scalar fres(const scalar ft, const scalar stoicRatio) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const ThermoType & mixture(const scalar, const scalar, const scalar) const
Definition egrMixture.C:60
const dimensionedScalar & stoicRatio() const
Definition egrMixture.H:122
void read(const dictionary &)
Read dictionary.
Definition egrMixture.C:91
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
Definition egrMixture.C:103
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for handling words, derived from Foam::string.
Definition word.H:66
PtrList< volScalarField > & Y
const dictionary & thermoDict
Definition EEqn.H:16
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const dimensionSet dimless
Dimensionless.
hashedWordList speciesTable
A table of species as a hashedWordList.
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
volScalarField & b