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