homogeneousMixture.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 "homogeneousMixture.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class ThermoType>
39 (
40  const dictionary& thermoDict,
41  const fvMesh& mesh,
42  const word& phaseName
43 )
44 :
46  (
47  thermoDict,
48  speciesTable(nSpecies_, specieNames_),
49  mesh,
50  phaseName
51  ),
52 
53  reactants_(thermoDict.subDict("reactants")),
54  products_(thermoDict.subDict("products")),
55  mixture_("mixture", reactants_),
56  b_(Y("b"))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class ThermoType>
64 (
65  const scalar b
66 ) const
67 {
68  if (b > 0.999)
69  {
70  return reactants_;
71  }
72  else if (b < 0.001)
73  {
74  return products_;
75  }
76  else
77  {
78  mixture_ = b*reactants_;
79  mixture_ += (1 - b)*products_;
80 
81  return mixture_;
82  }
83 }
84 
85 
86 template<class ThermoType>
88 {
89  reactants_ = ThermoType(thermoDict.subDict("reactants"));
90  products_ = ThermoType(thermoDict.subDict("products"));
91 }
92 
93 
94 template<class ThermoType>
96 (
97  const label speciei
98 ) const
99 {
100  if (speciei == 0)
101  {
102  return reactants_;
103  }
104  else if (speciei == 1)
105  {
106  return products_;
107  }
108  else
109  {
111  << "Unknown specie index " << speciei << ". Valid indices are 0..1"
112  << abort(FatalError);
113 
114  return reactants_;
115  }
116 }
117 
118 
119 // ************************************************************************* //
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:306
const thermoType & mixture(const scalar) const
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:1002
A class for handling words, derived from string.
Definition: word.H:59
Specialisation of the basicMixture for combustion.
const ThermoType & specieThermo(const label speciei) const
Return thermo based on index.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::homogeneousMixture.
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
void read(const dictionary &)
Read dictionary.
homogeneousMixture(const dictionary &, const fvMesh &, const word &)
Construct from dictionary, mesh and phase name.