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-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::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 (
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 (
87  const incompressiblePerfectGas<Specie>&,
88  const incompressiblePerfectGas<Specie>&
89 );
90 
91 template<class Specie>
92 Ostream& operator<<
93 (
94  Ostream&,
95  const incompressiblePerfectGas<Specie>&
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 dictionary
123 
124  //- Construct as named copy
126  (
127  const word& name,
129  );
130 
131  //- Construct and return a clone
133 
134  // Selector from dictionary
136  (
137  const dictionary& dict
138  );
139 
140 
141  // Member Functions
142 
143  //- Return the instantiated type name
144  static word typeName()
145  {
146  return
147  "incompressiblePerfectGas<"
148  + word(Specie::typeName_()) + '>';
149  }
150 
151 
152  // Fundamental properties
153 
154  //- Is the equation of state is incompressible i.e. rho != f(p)
155  static const bool incompressible = true;
156 
157  //- Is the equation of state is isochoric i.e. rho = const
158  static const bool isochoric = false;
159 
160  //- Return density [kg/m^3]
161  inline scalar rho(scalar p, scalar T) const;
162 
163  //- Return enthalpy contribution [J/kg]
164  inline scalar H(const scalar p, const scalar T) const;
165 
166  //- Return Cp contribution [J/(kg K]
167  inline scalar Cp(scalar p, scalar T) const;
168 
169  //- Return internal energy contribution [J/kg]
170  inline scalar E(const scalar p, const scalar T) const;
171 
172  //- Return Cv contribution [J/(kg K]
173  inline scalar Cv(scalar p, scalar T) const;
174 
175  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
176  inline scalar Sp(const scalar p, const scalar T) const;
177 
178  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
179  inline scalar Sv(const scalar p, const scalar T) const;
180 
181  //- Return compressibility [s^2/m^2]
182  inline scalar psi(scalar p, scalar T) const;
183 
184  //- Return compression factor []
185  inline scalar Z(scalar p, scalar T) const;
186 
187  //- Return (Cp - Cv) [J/(kg K]
188  inline scalar CpMCv(scalar p, scalar T) const;
189 
190  //- Return volumetric coefficient of thermal expansion [1/T]
191  inline scalar alphav(const scalar p, const scalar T) const;
192 
193 
194  // IO
195 
196  //- Write to Ostream
197  void write(Ostream& os) const;
198 
199 
200  // Member Operators
201 
202  inline void operator+=(const incompressiblePerfectGas&);
203  inline void operator*=(const scalar);
204 
205 
206  // Friend operators
207 
208  friend incompressiblePerfectGas operator+ <Specie>
209  (
211  const incompressiblePerfectGas&
212  );
213 
214  friend incompressiblePerfectGas operator* <Specie>
215  (
216  const scalar s,
217  const incompressiblePerfectGas&
218  );
219 
220  friend incompressiblePerfectGas operator== <Specie>
221  (
222  const incompressiblePerfectGas&,
223  const incompressiblePerfectGas&
224  );
225 
226 
227  // Ostream Operator
228 
229  friend Ostream& operator<< <Specie>
230  (
231  Ostream&,
232  const incompressiblePerfectGas&
233  );
234 };
235 
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 } // End namespace Foam
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
244 
245 #ifdef NoRepository
246  #include "incompressiblePerfectGas.C"
247 #endif
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
void write(Ostream &os) const
Write to Ostream.
dictionary dict
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Incompressible gas equation of state using a constant reference pressure in the perfect gas equation ...
void operator+=(const incompressiblePerfectGas &)
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].
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 Z(scalar p, scalar T) const
Return compression factor [].
A class for handling words, derived from string.
Definition: word.H:59
incompressiblePerfectGas(const Specie &sp, const scalar pRef)
Construct from components.
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
scalar pRef
Definition: createFields.H:6
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 H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [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
autoPtr< incompressiblePerfectGas > clone() const
Construct and return a clone.
static autoPtr< incompressiblePerfectGas > New(const dictionary &dict)
static word typeName()
Return the instantiated type name.
Namespace for OpenFOAM.