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