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