perfectGas.H
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-2021 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 Class
25  Foam::perfectGas
26 
27 Description
28  Perfect gas equation of state:
29 
30  \verbatim
31  rho = p/(R*T)
32  \endverbatim
33 
34 Usage
35  The gas constant R used by the perfect gas equation of state is obtained
36  directly from the specie or mixture molecular weight so there is no need to
37  provide an equation of \c equationOfState enry in the \c mixture
38  specification.
39 
40 SourceFiles
41  perfectGasI.H
42  perfectGas.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef perfectGas_H
47 #define perfectGas_H
48 
49 #include "autoPtr.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward declaration of friend functions and operators
57 
58 template<class Specie> class perfectGas;
59 
60 template<class Specie>
61 inline perfectGas<Specie> operator+
62 (
63  const perfectGas<Specie>&,
64  const perfectGas<Specie>&
65 );
66 
67 template<class Specie>
68 inline perfectGas<Specie> operator*
69 (
70  const scalar,
71  const perfectGas<Specie>&
72 );
73 
74 template<class Specie>
75 inline perfectGas<Specie> operator==
76 (
77  const perfectGas<Specie>&,
78  const perfectGas<Specie>&
79 );
80 
81 template<class Specie>
82 Ostream& operator<<
83 (
84  Ostream&,
85  const perfectGas<Specie>&
86 );
87 
88 
89 /*---------------------------------------------------------------------------*\
90  Class perfectGas Declaration
91 \*---------------------------------------------------------------------------*/
92 
93 template<class Specie>
94 class perfectGas
95 :
96  public Specie
97 {
98 
99 public:
100 
101  // Constructors
102 
103  //- Construct from components
104  inline perfectGas(const Specie& sp);
105 
106  //- Construct from dictionary
107  perfectGas(const dictionary& dict);
108 
109  //- Construct as named copy
110  inline perfectGas(const word& name, const perfectGas&);
111 
112  //- Construct and return a clone
113  inline autoPtr<perfectGas> clone() const;
114 
115  // Selector from dictionary
116  inline static autoPtr<perfectGas> New(const dictionary& dict);
117 
118 
119  // Member Functions
120 
121  //- Return the instantiated type name
122  static word typeName()
123  {
124  return "perfectGas<" + word(Specie::typeName_()) + '>';
125  }
126 
127 
128  // Fundamental properties
129 
130  //- Is the equation of state is incompressible i.e. rho != f(p)
131  static const bool incompressible = false;
132 
133  //- Is the equation of state is isochoric i.e. rho = const
134  static const bool isochoric = false;
135 
136  //- Return density [kg/m^3]
137  inline scalar rho(scalar p, scalar T) const;
138 
139  //- Return enthalpy contribution [J/kg]
140  inline scalar H(const scalar p, const scalar T) const;
141 
142  //- Return Cp contribution [J/(kg K]
143  inline scalar Cp(scalar p, scalar T) const;
144 
145  //- Return internal energy contribution [J/kg]
146  inline scalar E(const scalar p, const scalar T) const;
147 
148  //- Return Cv contribution [J/(kg K]
149  inline scalar Cv(scalar p, scalar T) const;
150 
151  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
152  inline scalar Sp(const scalar p, const scalar T) const;
153 
154  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
155  inline scalar Sv(const scalar p, const scalar T) const;
156 
157  //- Return compressibility [s^2/m^2]
158  inline scalar psi(scalar p, scalar T) const;
159 
160  //- Return compression factor []
161  inline scalar Z(scalar p, scalar T) const;
162 
163  //- Return (Cp - Cv) [J/(kg K]
164  inline scalar CpMCv(scalar p, scalar T) const;
165 
166  //- Return volumetric coefficient of thermal expansion [1/T]
167  inline scalar alphav(const scalar p, const scalar T) const;
168 
169 
170  // IO
171 
172  //- Write to Ostream
173  void write(Ostream& os) const;
174 
175 
176  // Member Operators
177 
178  inline void operator+=(const perfectGas&);
179  inline void operator*=(const scalar);
180 
181 
182  // Friend operators
183 
184  friend perfectGas operator+ <Specie>
185  (
186  const perfectGas&,
187  const perfectGas&
188  );
189 
190  friend perfectGas operator* <Specie>
191  (
192  const scalar s,
193  const perfectGas&
194  );
195 
196  friend perfectGas operator== <Specie>
197  (
198  const perfectGas&,
199  const perfectGas&
200  );
201 
202 
203  // Ostream Operator
204 
205  friend Ostream& operator<< <Specie>
206  (
207  Ostream&,
208  const perfectGas&
209  );
210 };
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #include "perfectGasI.H"
220 
221 #ifdef NoRepository
222  #include "perfectGas.C"
223 #endif
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 #endif
228 
229 // ************************************************************************* //
void operator*=(const scalar)
Definition: perfectGasI.H:158
static word typeName()
Return the instantiated type name.
Definition: perfectGas.H:121
dictionary dict
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: perfectGasI.H:135
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: perfectGasI.H:113
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: perfectGasI.H:106
void write(Ostream &os) const
Write to Ostream.
Definition: perfectGas.C:41
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: perfectGasI.H:99
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: perfectGasI.H:78
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: perfectGasI.H:85
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: perfectGasI.H:121
void operator+=(const perfectGas &)
Definition: perfectGasI.H:151
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: perfectGas.H:130
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))
autoPtr< perfectGas > clone() const
Construct and return a clone.
Definition: perfectGasI.H:52
A class for handling words, derived from string.
Definition: word.H:59
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: perfectGasI.H:92
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Perfect gas equation of state:
Definition: perfectGas.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static autoPtr< perfectGas > New(const dictionary &dict)
Definition: perfectGasI.H:60
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: perfectGas.H:133
perfectGas(const Specie &sp)
Construct from components.
Definition: perfectGasI.H:31
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectGasI.H:128
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectGasI.H:71
Namespace for OpenFOAM.
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: perfectGasI.H:142