SpecieMixture.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-2018 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 "SpecieMixture.H"
27 #include "fvMesh.H"
28 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class MixtureType>
34 (
35  const dictionary& thermoDict,
36  const fvMesh& mesh,
37  const word& phaseName
38 )
39 :
40  MixtureType
41  (
42  thermoDict,
43  mesh,
44  phaseName
45  )
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
51 template<class MixtureType>
53 (
54  const label speciei
55 ) const
56 {
57  return this->getLocalThermo(speciei).W();
58 }
59 
60 
61 template<class MixtureType>
63 (
64  const label speciei,
65  const scalar p,
66  const scalar T
67 ) const
68 {
69  return this->getLocalThermo(speciei).Cp(p, T);
70 }
71 
72 
73 template<class MixtureType>
75 (
76  const label speciei,
77  const scalar p,
78  const scalar T
79 ) const
80 {
81  return this->getLocalThermo(speciei).Cv(p, T);
82 }
83 
84 
85 template<class MixtureType>
87 (
88  const label speciei,
89  const scalar p,
90  const scalar T
91 ) const
92 {
93  return this->getLocalThermo(speciei).Ha(p, T);
94 }
95 
96 
97 template<class MixtureType>
99 (
100  const label speciei,
101  const scalar p,
102  const scalar T
103 ) const
104 {
105  return this->getLocalThermo(speciei).Hs(p, T);
106 }
107 
108 
109 template<class MixtureType>
111 (
112  const label speciei
113 ) const
114 {
115  return this->getLocalThermo(speciei).Hc();
116 }
117 
118 
119 template<class MixtureType>
121 (
122  const label speciei,
123  const scalar p,
124  const scalar T
125 ) const
126 {
127  return this->getLocalThermo(speciei).S(p, T);
128 }
129 
130 
131 template<class MixtureType>
133 (
134  const label speciei,
135  const scalar p,
136  const scalar T
137 ) const
138 {
139  return this->getLocalThermo(speciei).Es(p, T);
140 }
141 
142 
143 template<class MixtureType>
145 (
146  const label speciei,
147  const scalar p,
148  const scalar T
149 ) const
150 {
151  return this->getLocalThermo(speciei).G(p, T);
152 }
153 
154 
155 template<class MixtureType>
157 (
158  const label speciei,
159  const scalar p,
160  const scalar T
161 ) const
162 {
163  return this->getLocalThermo(speciei).A(p, T);
164 }
165 
166 
167 template<class MixtureType>
169 (
170  const label speciei,
171  const scalar p,
172  const scalar T
173 ) const
174 {
175  return this->getLocalThermo(speciei).mu(p, T);
176 }
177 
178 
179 template<class MixtureType>
181 (
182  const label speciei,
183  const scalar p,
184  const scalar T
185 ) const
186 {
187  return this->getLocalThermo(speciei).kappa(p, T);
188 }
189 
190 
191 template<class MixtureType>
193 (
194  const label speciei,
195  const scalar p,
196  const scalar T
197 ) const
198 {
199  return this->getLocalThermo(speciei).alphah(p, T);
200 }
201 
202 
203 template<class MixtureType>
205 (
206  const label speciei,
207  const scalar p,
208  const scalar T
209 ) const
210 {
211  return this->getLocalThermo(speciei).rho(p, T);
212 }
213 
214 
215 // ************************************************************************* //
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: SpecieMixture.C:99
virtual scalar Wi(const label speciei) const
Molecular weight of the given specie [kg/kmol].
Definition: SpecieMixture.C:53
virtual scalar S(const label speciei, const scalar p, const scalar T) const
Entropy [J/kg/K].
virtual scalar Hc(const label speciei) const
Chemical enthalpy [J/kg].
virtual scalar alphah(const label speciei, const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/m/s].
SpecieMixture(const dictionary &, const fvMesh &, const word &phaseName)
Construct from dictionary, mesh and phase name.
Definition: SpecieMixture.C:34
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< volScalarField > Cv() const =0
Heat capacity at constant volume [J/kg/K].
virtual scalar Ha(const label speciei, const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
Definition: SpecieMixture.C:87
virtual scalar A(const label speciei, const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
virtual scalar rho(const label speciei, const scalar p, const scalar T) const
Density [kg/m^3].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
virtual scalar Es(const label speciei, const scalar p, const scalar T) const
Sensible internal energy [J/kg].
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [W/m/K].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual scalar G(const label speciei, const scalar p, const scalar T) const
Gibbs free energy [J/kg].
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure for patch [J/kg/K].