perfectFluid.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2012-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::perfectFluid
28
29Group
30 grpSpecieEquationOfState
31
32Description
33 Perfect gas equation of state.
34
35SourceFiles
36 perfectFluidI.H
37 perfectFluid.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef perfectFluid_H
42#define perfectFluid_H
43
44#include "autoPtr.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51// Forward declaration of friend functions and operators
52
53template<class Specie> class perfectFluid;
54
55template<class Specie>
56inline perfectFluid<Specie> operator+
57(
60);
61
62template<class Specie>
63inline perfectFluid<Specie> operator*
64(
65 const scalar,
67);
68
69template<class Specie>
70inline perfectFluid<Specie> operator==
71(
74);
75
76template<class Specie>
77Ostream& operator<<
78(
79 Ostream&,
81);
82
83
84/*---------------------------------------------------------------------------*\
85 Class perfectFluid Declaration
86\*---------------------------------------------------------------------------*/
87
88template<class Specie>
89class perfectFluid
90:
91 public Specie
92{
93 // Private data
94
95 //- Fluid constant
96 scalar R_;
97
98 //- The reference density
99 scalar rho0_;
100
101
102public:
103
104 // Constructors
105
106 //- Construct from components
107 inline perfectFluid
108 (
109 const Specie& sp,
110 const scalar R,
111 const scalar rho0
112 );
113
114 //- Construct from dictionary
116
117 //- Construct as named copy
118 inline perfectFluid(const word& name, const perfectFluid&);
119
120 //- Construct and return a clone
121 inline autoPtr<perfectFluid> clone() const;
122
123 // Selector from dictionary
124 inline static autoPtr<perfectFluid> New(const dictionary& dict);
125
126
127 // Member functions
128
129 //- Return the instantiated type name
130 static word typeName()
131 {
132 return "perfectFluid<" + word(Specie::typeName_()) + '>';
133 }
134
135
136 // Fundamental properties
137
138 //- Is the equation of state is incompressible i.e. rho != f(p)
139 static const bool incompressible = false;
140
141 //- Is the equation of state is isochoric i.e. rho = const
142 static const bool isochoric = false;
143
144 //- Return fluid constant [J/(kg K)]
145 inline scalar R() const;
146
147 //- Return density [kg/m^3]
148 inline scalar rho(scalar p, scalar T) const;
149
150 //- Return enthalpy departure [J/kg]
151 inline scalar H(const scalar p, const scalar T) const;
152
153 //- Return Cp departure [J/(kg K]
154 inline scalar Cp(scalar p, scalar T) const;
155
156 //- Return internal energy departure [J/kg]
157 inline scalar E(const scalar p, const scalar T) const;
158
159 //- Return Cv departure [J/(kg K]
160 inline scalar Cv(scalar p, scalar T) const;
161
162 //- Return entropy [J/(kg K)]
163 inline scalar S(const scalar p, const scalar T) const;
164
165 //- Return compressibility rho/p [s^2/m^2]
166 inline scalar psi(scalar p, scalar T) const;
167
168 //- Return compression factor []
169 inline scalar Z(scalar p, scalar T) const;
170
171 //- Return (Cp - Cv) [J/(kg K]
172 inline scalar CpMCv(scalar p, scalar T) const;
173
174
175 // IO
176
177 //- Write to Ostream
178 void write(Ostream& os) const;
179
180
181 // Member operators
182
183 inline void operator+=(const perfectFluid&);
184 inline void operator*=(const scalar);
185
186
187 // Friend operators
189 friend perfectFluid operator+ <Specie>
190 (
191 const perfectFluid&,
192 const perfectFluid&
193 );
195 friend perfectFluid operator* <Specie>
196 (
197 const scalar s,
198 const perfectFluid&
199 );
201 friend perfectFluid operator== <Specie>
202 (
203 const perfectFluid&,
204 const perfectFluid&
205 );
206
207
208 // Ostream Operator
210 friend Ostream& operator<< <Specie>
211 (
212 Ostream&,
213 const perfectFluid&
214 );
215};
216
217
218// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219
220} // End namespace Foam
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
224#include "perfectFluidI.H"
225
226#ifdef NoRepository
227 #include "perfectFluid.C"
228#endif
229
230// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231
232#endif
233
234// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Perfect gas equation of state.
Definition: perfectFluid.H:91
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: perfectFluidI.H:98
void operator+=(const perfectFluid &)
static word typeName()
Return the instantiated type name.
Definition: perfectFluid.H:129
perfectFluid(const word &name, const perfectFluid &)
Construct as named copy.
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
autoPtr< perfectFluid > clone() const
Construct and return a clone.
Definition: perfectFluidI.H:64
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
static autoPtr< perfectFluid > New(const dictionary &dict)
Definition: perfectFluidI.H:73
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: perfectFluid.H:141
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: perfectFluid.H:138
scalar Z(scalar p, scalar T) const
Return compression factor [].
scalar R() const
Return fluid constant [J/(kg K)].
Definition: perfectFluidI.H:84
void operator*=(const scalar)
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & psi
const volScalarField & Cv
Definition: EEqn.H:8
const volScalarField & Cp
Definition: EEqn.H:7
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
runTime write()
scalar rho0
dictionary dict