perfectGasI.H
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-2016 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 "perfectGas.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
31 inline Foam::perfectGas<Specie>::perfectGas(const Specie& sp)
32 :
33  Specie(sp)
34 {}
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class Specie>
41 (
42  const word& name,
43  const perfectGas<Specie>& pg
44 )
45 :
46  Specie(name, pg)
47 {}
48 
49 
50 template<class Specie>
53 {
55 }
56 
57 
58 template<class Specie>
61 {
63 }
64 
65 
66 template<class Specie>
68 (
69  const dictionary& dict
70 )
71 {
73 }
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
78 template<class Specie>
79 inline Foam::scalar Foam::perfectGas<Specie>::rho(scalar p, scalar T) const
80 {
81  return p/(this->R()*T);
82 }
83 
84 
85 template<class Specie>
86 inline Foam::scalar Foam::perfectGas<Specie>::h(scalar p, scalar T) const
87 {
88  return 0;
89 }
90 
91 
92 template<class Specie>
93 inline Foam::scalar Foam::perfectGas<Specie>::cp(scalar p, scalar T) const
94 {
95  return 0;
96 }
97 
98 
99 template<class Specie>
100 inline Foam::scalar Foam::perfectGas<Specie>::s(scalar p, scalar T) const
101 {
102  return -RR*log(p/Pstd);
103 }
104 
105 
106 template<class Specie>
107 inline Foam::scalar Foam::perfectGas<Specie>::psi(scalar p, scalar T) const
108 {
109  return 1.0/(this->R()*T);
110 }
111 
112 
113 template<class Specie>
114 inline Foam::scalar Foam::perfectGas<Specie>::Z(scalar p, scalar T) const
115 {
116  return 1;
117 }
118 
119 
120 template<class Specie>
121 inline Foam::scalar Foam::perfectGas<Specie>::cpMcv(scalar p, scalar T) const
122 {
123  return RR;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
128 
129 template<class Specie>
131 {
132  Specie::operator+=(pg);
133 }
134 
135 
136 template<class Specie>
138 {
139  Specie::operator-=(pg);
140 }
141 
142 
143 template<class Specie>
144 inline void Foam::perfectGas<Specie>::operator*=(const scalar s)
145 {
146  Specie::operator*=(s);
147 }
148 
149 
150 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
151 
152 template<class Specie>
153 inline Foam::perfectGas<Specie> Foam::operator+
154 (
155  const perfectGas<Specie>& pg1,
156  const perfectGas<Specie>& pg2
157 )
158 {
159  return perfectGas<Specie>
160  (
161  static_cast<const Specie&>(pg1)
162  + static_cast<const Specie&>(pg2)
163  );
164 }
165 
166 
167 template<class Specie>
168 inline Foam::perfectGas<Specie> Foam::operator-
169 (
170  const perfectGas<Specie>& pg1,
171  const perfectGas<Specie>& pg2
172 )
173 {
174  return perfectGas<Specie>
175  (
176  static_cast<const Specie&>(pg1)
177  - static_cast<const Specie&>(pg2)
178  );
179 }
180 
181 
182 template<class Specie>
183 inline Foam::perfectGas<Specie> Foam::operator*
184 (
185  const scalar s,
186  const perfectGas<Specie>& pg
187 )
188 {
189  return perfectGas<Specie>(s*static_cast<const Specie&>(pg));
190 }
191 
192 
193 template<class Specie>
194 inline Foam::perfectGas<Specie> Foam::operator==
195 (
196  const perfectGas<Specie>& pg1,
197  const perfectGas<Specie>& pg2
198 )
199 {
200  return pg2 - pg1;
201 }
202 
203 
204 // ************************************************************************* //
void operator*=(const scalar)
Definition: perfectGasI.H:144
static autoPtr< perfectGas > New(Istream &is)
Definition: perfectGasI.H:60
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectGasI.H:79
void operator+=(const perfectGas &)
Definition: perfectGasI.H:130
A class for handling words, derived from string.
Definition: word.H:59
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Definition: perfectGasI.H:100
const dimensionedScalar Pstd
Standard pressure.
void operator-=(const perfectGas &)
Definition: perfectGasI.H:137
Perfect gas equation of state.
Definition: perfectGas.H:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< perfectGas > clone() const
Construct and return a clone.
Definition: perfectGasI.H:52
scalar h(const scalar p, const scalar T) const
Return enthalpy departure [J/kmol].
Definition: perfectGasI.H:86
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Definition: perfectGasI.H:107
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
Definition: perfectGasI.H:121
#define R(A, B, C, D, E, F, K, M)
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectGasI.H:114
perfectGas(const Specie &sp)
Construct from components.
Definition: perfectGasI.H:31
const scalar RR
Universal gas constant (default in [J/(kmol K)])
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
volScalarField & p
scalar cp(scalar p, scalar T) const
Return cp departure [J/(kmol K].
Definition: perfectGasI.H:93