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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-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::adiabaticPerfectFluid
28
29Group
30 grpSpecieEquationOfState
31
32Description
33 Adiabatic perfect fluid equation of state.
34
35SourceFiles
36 adiabaticPerfectFluidI.H
37 adiabaticPerfectFluid.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef adiabaticPerfectFluid_H
42#define adiabaticPerfectFluid_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 adiabaticPerfectFluid;
54
55template<class Specie>
56inline adiabaticPerfectFluid<Specie> operator+
57(
60);
61
62template<class Specie>
63inline adiabaticPerfectFluid<Specie> operator*
64(
65 const scalar,
67);
68
69template<class Specie>
70inline adiabaticPerfectFluid<Specie> operator==
71(
74);
75
76template<class Specie>
77Ostream& operator<<
78(
79 Ostream&,
81);
82
83
84/*---------------------------------------------------------------------------*\
85 Class adiabaticPerfectFluid Declaration
86\*---------------------------------------------------------------------------*/
87
88template<class Specie>
90:
91 public Specie
92{
93 // Private data
94
95 //- Reference pressure
96 scalar p0_;
97
98 //- Reference density
99 scalar rho0_;
100
101 //- The isentropic exponent
102 scalar gamma_;
103
104 //- Pressure offset for a stiffened gas
105 scalar B_;
106
107
108public:
109
110 // Constructors
111
112 //- Construct from components
114 (
115 const Specie& sp,
116 const scalar p0,
117 const scalar rho0,
118 const scalar gamma,
119 const scalar B
120 );
121
122 //- Construct from dictionary
124
125 //- Construct as named copy
127 (
128 const word& name,
130 );
131
132 //- Construct and return a clone
134
135 // Selector from dictionary
137 (
138 const dictionary& dict
139 );
140
141
142 // Member functions
143
144 //- Return the instantiated type name
145 static word typeName()
146 {
147 return "adiabaticPerfectFluid<" + word(Specie::typeName_()) + '>';
148 }
149
150
151 // Fundamental properties
152
153 //- Is the equation of state is incompressible i.e. rho != f(p)
154 static const bool incompressible = false;
155
156 //- Is the equation of state is isochoric i.e. rho = const
157 static const bool isochoric = false;
158
159 //- Return density [kg/m^3]
160 inline scalar rho(scalar p, scalar T) const;
161
162 //- Return enthalpy departure [J/kg]
163 inline scalar H(const scalar p, const scalar T) const;
164
165 //- Return Cp departure [J/(kg K]
166 inline scalar Cp(scalar p, scalar T) const;
167
168 //- Return internal energy departure [J/kg]
169 inline scalar E(const scalar p, const scalar T) const;
170
171 //- Return Cv departure [J/(kg K]
172 inline scalar Cv(scalar p, scalar T) const;
173
174 //- Return entropy [J/(kg K)]
175 inline scalar S(const scalar p, const scalar T) const;
176
177 //- Return compressibility rho/p [s^2/m^2]
178 inline scalar psi(scalar p, scalar T) const;
179
180 //- Return compression factor []
181 inline scalar Z(scalar p, scalar T) const;
182
183 //- Return (Cp - Cv) [J/(kg K]
184 inline scalar CpMCv(scalar p, scalar T) const;
185
186
187 // IO
188
189 //- Write to Ostream
190 void write(Ostream& os) const;
191
192
193 // Member operators
194
195 inline void operator+=(const adiabaticPerfectFluid&);
196 inline void operator*=(const scalar);
197
198
199 // Friend operators
201 friend adiabaticPerfectFluid operator+ <Specie>
202 (
205 );
207 friend adiabaticPerfectFluid operator* <Specie>
208 (
209 const scalar s,
211 );
213 friend adiabaticPerfectFluid operator== <Specie>
214 (
217 );
218
219
220 // Ostream Operator
222 friend Ostream& operator<< <Specie>
223 (
224 Ostream&,
226 );
227};
228
229
230// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231
232} // End namespace Foam
233
234// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235
237
238#ifdef NoRepository
239 #include "adiabaticPerfectFluid.C"
240#endif
241
242// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243
244#endif
245
246// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Adiabatic perfect fluid equation of state.
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].
static word typeName()
Return the instantiated type name.
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
static autoPtr< adiabaticPerfectFluid > New(const dictionary &dict)
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
adiabaticPerfectFluid(const word &name, const adiabaticPerfectFluid &)
Construct as named copy.
void operator+=(const adiabaticPerfectFluid &)
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
scalar Z(scalar p, scalar T) const
Return compression factor [].
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
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & psi
const volScalarField & T
const volScalarField & Cv
Definition: EEqn.H:8
const scalar gamma
Definition: EEqn.H:9
const volScalarField & Cp
Definition: EEqn.H:7
const volScalarField & p0
Definition: EEqn.H:36
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