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-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 "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>::s(scalar p, scalar T) const
87 {
88  return -RR*log(p/Pstd);
89 }
90 
91 
92 template<class Specie>
93 inline Foam::scalar Foam::perfectGas<Specie>::psi(scalar p, scalar T) const
94 {
95  return 1.0/(this->R()*T);
96 }
97 
98 
99 template<class Specie>
100 inline Foam::scalar Foam::perfectGas<Specie>::Z(scalar p, scalar T) const
101 {
102  return 1;
103 }
104 
105 
106 template<class Specie>
107 inline Foam::scalar Foam::perfectGas<Specie>::cpMcv(scalar p, scalar T) const
108 {
109  return RR;
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
114 
115 template<class Specie>
117 {
118  Specie::operator+=(pg);
119 }
120 
121 
122 template<class Specie>
124 {
125  Specie::operator-=(pg);
126 }
127 
128 
129 template<class Specie>
130 inline void Foam::perfectGas<Specie>::operator*=(const scalar s)
131 {
132  Specie::operator*=(s);
133 }
134 
135 
136 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
137 
138 template<class Specie>
139 inline Foam::perfectGas<Specie> Foam::operator+
140 (
141  const perfectGas<Specie>& pg1,
142  const perfectGas<Specie>& pg2
143 )
144 {
145  return perfectGas<Specie>
146  (
147  static_cast<const Specie&>(pg1)
148  + static_cast<const Specie&>(pg2)
149  );
150 }
151 
152 
153 template<class Specie>
154 inline Foam::perfectGas<Specie> Foam::operator-
155 (
156  const perfectGas<Specie>& pg1,
157  const perfectGas<Specie>& pg2
158 )
159 {
160  return perfectGas<Specie>
161  (
162  static_cast<const Specie&>(pg1)
163  - static_cast<const Specie&>(pg2)
164  );
165 }
166 
167 
168 template<class Specie>
169 inline Foam::perfectGas<Specie> Foam::operator*
170 (
171  const scalar s,
172  const perfectGas<Specie>& pg
173 )
174 {
175  return perfectGas<Specie>(s*static_cast<const Specie&>(pg));
176 }
177 
178 
179 template<class Specie>
180 inline Foam::perfectGas<Specie> Foam::operator==
181 (
182  const perfectGas<Specie>& pg1,
183  const perfectGas<Specie>& pg2
184 )
185 {
186  return pg2 - pg1;
187 }
188 
189 
190 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Definition: thermodynamicConstants.C:44
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constant::standard::Pstd
const dimensionedScalar Pstd
Standard pressure.
Definition: thermodynamicConstants.C:46
Foam::perfectGas::Z
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectGasI.H:100
Foam::perfectGas::perfectGas
perfectGas(const Specie &sp)
Construct from components.
Definition: perfectGasI.H:31
Foam::perfectGas::New
static autoPtr< perfectGas > New(Istream &is)
Definition: perfectGasI.H:60
R
#define R(A, B, C, D, E, F, K, M)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::perfectGas::s
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Definition: perfectGasI.H:86
Foam::perfectGas
Perfect gas equation of state.
Definition: perfectGas.H:47
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::perfectGas::clone
autoPtr< perfectGas > clone() const
Construct and return a clone.
Definition: perfectGasI.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
perfectGas.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
Foam::perfectGas::operator-=
void operator-=(const perfectGas &)
Definition: perfectGasI.H:123
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::perfectGas::operator+=
void operator+=(const perfectGas &)
Definition: perfectGasI.H:116
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::perfectGas::cpMcv
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
Definition: perfectGasI.H:107
Foam::perfectGas::psi
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Definition: perfectGasI.H:93
Foam::perfectGas::rho
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectGasI.H:79
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::perfectGas::operator*=
void operator*=(const scalar)
Definition: perfectGasI.H:130