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