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