Loading...
Searching...
No Matches
inhomogeneousMixture.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\*---------------------------------------------------------------------------*/
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class ThermoType>
34Foam::inhomogeneousMixture<ThermoType>::inhomogeneousMixture
35(
37 const fvMesh& mesh,
38 const word& phaseName
39)
40:
42 (
44 speciesTable({"ft", "b"}),
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{}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
64template<class ThermoType>
66(
67 const scalar ft,
68 const scalar b
69) const
70{
71 if (ft < 0.0001)
72 {
73 return oxidant_;
74 }
75 else
76 {
77 scalar fu = b*ft + (1.0 - b)*fres(ft, stoicRatio().value());
78 scalar ox = 1 - ft - (ft - fu)*stoicRatio().value();
79 scalar pr = 1 - fu - ox;
80
81 mixture_ = fu*fuel_;
82 mixture_ += ox*oxidant_;
83 mixture_ += pr*products_;
85 return mixture_;
86 }
87}
88
89
90template<class ThermoType>
92{
93 thermoDict.readEntry("stoichiometricAirFuelMassRatio", stoicRatio_);
94
95 fuel_ = ThermoType(thermoDict.subDict("fuel"));
96 oxidant_ = ThermoType(thermoDict.subDict("oxidant"));
97 products_ = ThermoType(thermoDict.subDict("burntProducts"));
98}
99
100
101template<class ThermoType>
103(
104 const label speciei
105) const
106{
107 if (speciei == 0)
108 {
109 return fuel_;
110 }
111 else if (speciei == 1)
112 {
113 return oxidant_;
114 }
115 else if (speciei == 2)
116 {
117 return products_;
118 }
119 else
120 {
122 << "Unknown specie index " << speciei << ". Valid indices are 0..2"
123 << abort(FatalError);
124
125 return fuel_;
126 }
127}
128
129
130// ************************************************************************* //
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const ThermoType & mixture(const scalar, const scalar) const
const dimensionedScalar & stoicRatio() const
void read(const dictionary &)
Read dictionary.
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
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