SpecieMixture.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 | Copyright (C) 2011-2015 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).nMoles();
58 }
59 
60 
61 template<class MixtureType>
63 (
64  const label speciei
65 ) const
66 {
67  return this->getLocalThermo(speciei).W();
68 }
69 
70 
71 template<class MixtureType>
73 (
74  const label speciei,
75  const scalar p,
76  const scalar T
77 ) const
78 {
79  return this->getLocalThermo(speciei).Cp(p, T);
80 }
81 
82 
83 template<class MixtureType>
85 (
86  const label speciei,
87  const scalar p,
88  const scalar T
89 ) const
90 {
91  return this->getLocalThermo(speciei).Cv(p, T);
92 }
93 
94 
95 template<class MixtureType>
97 (
98  const label speciei,
99  const scalar p,
100  const scalar T
101 ) const
102 {
103  return this->getLocalThermo(speciei).Ha(p, T);
104 }
105 
106 
107 template<class MixtureType>
109 (
110  const label speciei,
111  const scalar p,
112  const scalar T
113 ) const
114 {
115  return this->getLocalThermo(speciei).Hs(p, T);
116 }
117 
118 
119 template<class MixtureType>
121 (
122  const label speciei
123 ) const
124 {
125  return this->getLocalThermo(speciei).Hc();
126 }
127 
128 
129 template<class MixtureType>
131 (
132  const label speciei,
133  const scalar p,
134  const scalar T
135 ) const
136 {
137  return this->getLocalThermo(speciei).S(p, T);
138 }
139 
140 
141 template<class MixtureType>
143 (
144  const label speciei,
145  const scalar p,
146  const scalar T
147 ) const
148 {
149  return this->getLocalThermo(speciei).Es(p, T);
150 }
151 
152 
153 template<class MixtureType>
155 (
156  const label speciei,
157  const scalar p,
158  const scalar T
159 ) const
160 {
161  return this->getLocalThermo(speciei).G(p, T);
162 }
163 
164 
165 template<class MixtureType>
167 (
168  const label speciei,
169  const scalar p,
170  const scalar T
171 ) const
172 {
173  return this->getLocalThermo(speciei).A(p, T);
174 }
175 
176 
177 template<class MixtureType>
179 (
180  const label speciei,
181  const scalar p,
182  const scalar T
183 ) const
184 {
185  return this->getLocalThermo(speciei).mu(p, T);
186 }
187 
188 
189 template<class MixtureType>
191 (
192  const label speciei,
193  const scalar p,
194  const scalar T
195 ) const
196 {
197  return this->getLocalThermo(speciei).kappa(p, T);
198 }
199 
200 
201 template<class MixtureType>
203 (
204  const label speciei,
205  const scalar p,
206  const scalar T
207 ) const
208 {
209  return this->getLocalThermo(speciei).alphah(p, T);
210 }
211 
212 
213 template<class MixtureType>
215 (
216  const label speciei,
217  const scalar p,
218  const scalar T
219 ) const
220 {
221  return this->getLocalThermo(speciei).rho(p, T);
222 }
223 
224 
225 // ************************************************************************* //
virtual scalar A(const label speciei, const scalar p, const scalar T) const
Helmholtz free energy [J/kg].
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
virtual scalar W(const label speciei) const
Molecular weight of the given specie [kg/kmol].
Definition: SpecieMixture.C:63
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual scalar alphah(const label speciei, const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/m/s].
virtual scalar rho(const label speciei, const scalar p, const scalar T) const
Density [kg/m3].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const
Dynamic viscosity [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 scalar G(const label speciei, const scalar p, const scalar T) const
Gibbs free energy [J/kg].
virtual scalar nMoles(const label speciei) const
Number of moles of the given specie [].
Definition: SpecieMixture.C:53
virtual tmp< volScalarField > Cv() const =0
Heat capacity at constant volume [J/kg/K].
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
virtual scalar Hc(const label speciei) const
Chemical enthalpy [J/kg].
virtual scalar Ha(const label speciei, const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
Definition: SpecieMixture.C:97
virtual scalar S(const label speciei, const scalar p, const scalar T) const
Entropy [J/(kg K)].
virtual scalar Es(const label speciei, const scalar p, const scalar T) const
Sensible internal energy [J/kg].
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure for patch [J/kg/K].