homogeneousMixture.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\*---------------------------------------------------------------------------*/
27
28#include "homogeneousMixture.H"
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class ThermoType>
35(
37 const fvMesh& mesh,
38 const word& phaseName
39)
40:
42 (
44 speciesTable({"b"}),
45 mesh,
46 phaseName
47 ),
48
49 reactants_(thermoDict.subDict("reactants")),
50 products_(thermoDict.subDict("products")),
51 mixture_("mixture", reactants_),
52 b_(Y("b"))
53{}
54
55
56// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57
58template<class ThermoType>
60(
61 const scalar b
62) const
63{
64 if (b > 0.999)
65 {
66 return reactants_;
67 }
68 else if (b < 0.001)
69 {
70 return products_;
71 }
72 else
73 {
74 mixture_ = b*reactants_;
75 mixture_ += (1 - b)*products_;
76
77 return mixture_;
78 }
79}
80
81
82template<class ThermoType>
84{
85 reactants_ = ThermoType(thermoDict.subDict("reactants"));
86 products_ = ThermoType(thermoDict.subDict("products"));
87}
88
89
90template<class ThermoType>
92(
93 const label speciei
94) const
95{
96 if (speciei == 0)
97 {
98 return reactants_;
99 }
100 else if (speciei == 1)
101 {
102 return products_;
103 }
104 else
105 {
107 << "Unknown specie index " << speciei << ". Valid indices are 0..1"
108 << abort(FatalError);
109
110 return reactants_;
111 }
112}
113
114
115// ************************************************************************* //
virtual bool read()
Re-read model coefficients if they have changed.
Specialization of the basicSpecieMixture for combustion.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
The homogeneous mixture contains species ("b").
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
A class for handling words, derived from Foam::string.
Definition: word.H:68
PtrList< volScalarField > & Y
const dictionary & thermoDict
Definition: EEqn.H:16
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
volScalarField & b
Definition: createFields.H:27