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-2018 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 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 scalar,
61  const adiabaticPerfectFluid<Specie>&
62 );
63 
64 template<class Specie>
65 inline adiabaticPerfectFluid<Specie> operator==
66 (
67  const adiabaticPerfectFluid<Specie>&,
68  const adiabaticPerfectFluid<Specie>&
69 );
70 
71 template<class Specie>
72 Ostream& operator<<
73 (
74  Ostream&,
75  const adiabaticPerfectFluid<Specie>&
76 );
77 
78 
79 /*---------------------------------------------------------------------------*\
80  Class adiabaticPerfectFluid Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 template<class Specie>
85 :
86  public Specie
87 {
88  // Private data
89 
90  //- Reference pressure
91  scalar p0_;
92 
93  //- Reference density
94  scalar rho0_;
95 
96  //- The isentropic exponent
97  scalar gamma_;
98 
99  //- Pressure offset for a stiffened gas
100  scalar B_;
101 
102 
103 public:
104 
105  // Constructors
106 
107  //- Construct from components
108  inline adiabaticPerfectFluid
109  (
110  const Specie& sp,
111  const scalar p0,
112  const scalar rho0,
113  const scalar gamma,
114  const scalar B
115  );
116 
117  //- Construct from dictionary
119 
120  //- Construct as named copy
121  inline adiabaticPerfectFluid
122  (
123  const word& name,
124  const adiabaticPerfectFluid&
125  );
126 
127  //- Construct and return a clone
128  inline autoPtr<adiabaticPerfectFluid> clone() const;
129 
130  // Selector from dictionary
131  inline static autoPtr<adiabaticPerfectFluid> New
132  (
133  const dictionary& dict
134  );
135 
136 
137  // Member functions
138 
139  //- Return the instantiated type name
140  static word typeName()
141  {
142  return "adiabaticPerfectFluid<" + word(Specie::typeName_()) + '>';
143  }
144 
145 
146  // Fundamental properties
147 
148  //- Is the equation of state is incompressible i.e. rho != f(p)
149  static const bool incompressible = false;
150 
151  //- Is the equation of state is isochoric i.e. rho = const
152  static const bool isochoric = false;
153 
154  //- Return density [kg/m^3]
155  inline scalar rho(scalar p, scalar T) const;
156 
157  //- Return enthalpy departure [J/kg]
158  inline scalar H(const scalar p, const scalar T) const;
159 
160  //- Return Cp departure [J/(kg K]
161  inline scalar Cp(scalar p, scalar T) const;
162 
163  //- Return entropy [J/(kg K)]
164  inline scalar S(const scalar p, const scalar T) const;
165 
166  //- Return compressibility rho/p [s^2/m^2]
167  inline scalar psi(scalar p, scalar T) const;
168 
169  //- Return compression factor []
170  inline scalar Z(scalar p, scalar T) const;
171 
172  //- Return (Cp - Cv) [J/(kg K]
173  inline scalar CpMCv(scalar p, scalar T) const;
174 
175 
176  // IO
177 
178  //- Write to Ostream
179  void write(Ostream& os) const;
180 
181 
182  // Member operators
183 
184  inline void operator+=(const adiabaticPerfectFluid&);
185  inline void operator*=(const scalar);
186 
187 
188  // Friend operators
189 
190  friend adiabaticPerfectFluid operator+ <Specie>
191  (
192  const adiabaticPerfectFluid&,
193  const adiabaticPerfectFluid&
194  );
195 
196  friend adiabaticPerfectFluid operator* <Specie>
197  (
198  const scalar s,
199  const adiabaticPerfectFluid&
200  );
201 
202  friend adiabaticPerfectFluid operator== <Specie>
203  (
204  const adiabaticPerfectFluid&,
205  const adiabaticPerfectFluid&
206  );
207 
208 
209  // Ostream Operator
210 
211  friend Ostream& operator<< <Specie>
212  (
213  Ostream&,
214  const adiabaticPerfectFluid&
215  );
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace Foam
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #include "adiabaticPerfectFluidI.H"
226 
227 #ifdef NoRepository
228  #include "adiabaticPerfectFluid.C"
229 #endif
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 #endif
234 
235 // ************************************************************************* //
scalar Z(scalar p, scalar T) const
Return compression factor [].
dictionary dict
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:137
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 rho/p [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 &)
scalar S(const scalar p, const scalar T) const
Return entropy [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 Cp(scalar p, scalar T) const
Return Cp departure [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 H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
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].
Adiabatic perfect fluid equation of state.
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.