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-2020 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.
29 
30  Coefficient mixing is very inaccurate and not supported,
31  so this equation of state is not applicable to mixtures.
32 
33 SourceFiles
34  adiabaticPerfectFluidI.H
35  adiabaticPerfectFluid.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef adiabaticPerfectFluid_H
40 #define adiabaticPerfectFluid_H
41 
42 #include "autoPtr.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // Forward declaration of friend functions and operators
50 
51 template<class Specie> class adiabaticPerfectFluid;
52 
53 template<class Specie>
54 inline adiabaticPerfectFluid<Specie> operator+
55 (
57  const adiabaticPerfectFluid<Specie>&
58 );
59 
60 template<class Specie>
61 inline adiabaticPerfectFluid<Specie> operator*
62 (
63  const scalar,
64  const adiabaticPerfectFluid<Specie>&
65 );
66 
67 template<class Specie>
68 inline adiabaticPerfectFluid<Specie> operator==
69 (
70  const adiabaticPerfectFluid<Specie>&,
71  const adiabaticPerfectFluid<Specie>&
72 );
73 
74 template<class Specie>
75 Ostream& operator<<
76 (
77  Ostream&,
78  const adiabaticPerfectFluid<Specie>&
79 );
80 
81 
82 /*---------------------------------------------------------------------------*\
83  Class adiabaticPerfectFluid Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 template<class Specie>
88 :
89  public Specie
90 {
91  // Private Data
92 
93  //- Reference pressure
94  scalar p0_;
95 
96  //- Reference density
97  scalar rho0_;
98 
99  //- The isentropic exponent
100  scalar gamma_;
101 
102  //- Pressure offset for a stiffened gas
103  scalar B_;
104 
105 
106 public:
107 
108  // Constructors
109 
110  //- Construct from components
111  inline adiabaticPerfectFluid
112  (
113  const Specie& sp,
114  const scalar p0,
115  const scalar rho0,
116  const scalar gamma,
117  const scalar B
118  );
119 
120  //- Construct from dictionary
122 
123  //- Construct as named copy
124  inline adiabaticPerfectFluid
125  (
126  const word& name,
127  const adiabaticPerfectFluid&
128  );
129 
130  //- Construct and return a clone
131  inline autoPtr<adiabaticPerfectFluid> clone() const;
132 
133  // Selector from dictionary
134  inline static autoPtr<adiabaticPerfectFluid> New
135  (
136  const dictionary& dict
137  );
138 
139 
140  // Member Functions
141 
142  //- Return the instantiated type name
143  static word typeName()
144  {
145  return "adiabaticPerfectFluid<" + word(Specie::typeName_()) + '>';
146  }
147 
148 
149  // Fundamental properties
150 
151  //- Is the equation of state is incompressible i.e. rho != f(p)
152  static const bool incompressible = false;
153 
154  //- Is the equation of state is isochoric i.e. rho = const
155  static const bool isochoric = false;
156 
157  //- Return density [kg/m^3]
158  inline scalar rho(scalar p, scalar T) const;
159 
160  //- Return enthalpy contribution [J/kg]
161  inline scalar H(const scalar p, const scalar T) const;
162 
163  //- Return Cp contribution [J/(kg K]
164  inline scalar Cp(scalar p, scalar T) const;
165 
166  //- Return internal energy contribution [J/kg]
167  inline scalar E(const scalar p, const scalar T) const;
168 
169  //- Return Cv contribution [J/(kg K]
170  inline scalar Cv(scalar p, scalar T) const;
171 
172  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
173  inline scalar Sp(const scalar p, const scalar T) const;
174 
175  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
176  inline scalar Sv(const scalar p, const scalar T) const;
177 
178  //- Return compressibility [s^2/m^2]
179  inline scalar psi(scalar p, scalar T) const;
180 
181  //- Return compression factor []
182  inline scalar Z(scalar p, scalar T) const;
183 
184  //- Return (Cp - Cv) [J/(kg K]
185  inline scalar CpMCv(scalar p, 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 adiabaticPerfectFluid&);
197  inline void operator*=(const scalar);
198 
199 
200  // Friend operators
201 
202  friend adiabaticPerfectFluid operator+ <Specie>
203  (
204  const adiabaticPerfectFluid&,
205  const adiabaticPerfectFluid&
206  );
207 
208  friend adiabaticPerfectFluid operator* <Specie>
209  (
210  const scalar s,
211  const adiabaticPerfectFluid&
212  );
213 
214  friend adiabaticPerfectFluid operator== <Specie>
215  (
216  const adiabaticPerfectFluid&,
217  const adiabaticPerfectFluid&
218  );
219 
220 
221  // Ostream Operator
222 
223  friend Ostream& operator<< <Specie>
224  (
225  Ostream&,
226  const adiabaticPerfectFluid&
227  );
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #include "adiabaticPerfectFluidI.H"
238 
239 #ifdef NoRepository
240  #include "adiabaticPerfectFluid.C"
241 #endif
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #endif
246 
247 // ************************************************************************* //
scalar Z(scalar p, scalar T) const
Return compression factor [].
dictionary dict
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
static autoPtr< adiabaticPerfectFluid > New(const dictionary &dict)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
adiabaticPerfectFluid(const Specie &sp, const scalar p0, const scalar rho0, const scalar gamma, const scalar B)
Construct from components.
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
static word typeName()
Return the instantiated type name.
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
scalar rho0
void operator+=(const adiabaticPerfectFluid &)
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))
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
A class for handling words, derived from string.
Definition: word.H:59
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
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
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Adiabatic perfect fluid equation of state.
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
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
Namespace for OpenFOAM.
void write(Ostream &os) const
Write to Ostream.