perfectGas.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) 2011-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::perfectGas
28
29Group
30 grpSpecieEquationOfState
31
32Description
33 Perfect gas equation of state.
34
35SourceFiles
36 perfectGasI.H
37 perfectGas.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef perfectGas_H
42#define perfectGas_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 perfectGas;
54
55template<class Specie>
56inline perfectGas<Specie> operator+
57(
58 const perfectGas<Specie>&,
60);
61
62template<class Specie>
63inline perfectGas<Specie> operator*
64(
65 const scalar,
67);
68
69template<class Specie>
70inline perfectGas<Specie> operator==
71(
72 const perfectGas<Specie>&,
74);
75
76template<class Specie>
77Ostream& operator<<
78(
79 Ostream&,
81);
82
83
84/*---------------------------------------------------------------------------*\
85 Class perfectGas Declaration
86\*---------------------------------------------------------------------------*/
87
88template<class Specie>
89class perfectGas
90:
91 public Specie
92{
93
94public:
95
96 // Constructors
97
98 //- Construct from components
99 inline perfectGas(const Specie& sp);
100
101 //- Construct from dictionary
102 perfectGas(const dictionary& dict);
103
104 //- Construct as named copy
105 inline perfectGas(const word& name, const perfectGas&);
106
107 //- Construct and return a clone
108 inline autoPtr<perfectGas> clone() const;
109
110 // Selector from dictionary
111 inline static autoPtr<perfectGas> New(const dictionary& dict);
112
113
114 // Member functions
115
116 //- Return the instantiated type name
117 static word typeName()
118 {
119 return "perfectGas<" + word(Specie::typeName_()) + '>';
120 }
121
122
123 // Fundamental properties
124
125 //- Is the equation of state is incompressible i.e. rho != f(p)
126 static const bool incompressible = false;
127
128 //- Is the equation of state is isochoric i.e. rho = const
129 static const bool isochoric = false;
130
131 //- Return density [kg/m^3]
132 inline scalar rho(scalar p, scalar T) const;
133
134 //- Return enthalpy departure [J/kg]
135 inline scalar H(const scalar p, const scalar T) const;
136
137 //- Return Cp departure [J/(kg K]
138 inline scalar Cp(scalar p, scalar T) const;
139
140 //- Return internal energy departure [J/kg]
141 inline scalar E(const scalar p, const scalar T) const;
142
143 //- Return Cv departure [J/(kg K]
144 inline scalar Cv(scalar p, scalar T) const;
145
146 //- Return entropy [J/(kg K)]
147 inline scalar S(const scalar p, const scalar T) const;
148
149 //- Return compressibility rho/p [s^2/m^2]
150 inline scalar psi(scalar p, scalar T) const;
151
152 //- Return compression factor []
153 inline scalar Z(scalar p, scalar T) const;
154
155 //- Return (Cp - Cv) [J/(kg K]
156 inline scalar CpMCv(scalar p, scalar T) const;
157
158
159 // IO
160
161 //- Write to Ostream
162 void write(Ostream& os) const;
163
164
165 // Member operators
166
167 inline void operator+=(const perfectGas&);
168 inline void operator*=(const scalar);
169
170
171 // Friend operators
173 friend perfectGas operator+ <Specie>
174 (
175 const perfectGas&,
176 const perfectGas&
177 );
179 friend perfectGas operator* <Specie>
180 (
181 const scalar s,
182 const perfectGas&
183 );
185 friend perfectGas operator== <Specie>
186 (
187 const perfectGas&,
188 const perfectGas&
189 );
190
191
192 // Ostream Operator
194 friend Ostream& operator<< <Specie>
195 (
196 Ostream&,
197 const perfectGas&
198 );
199};
200
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203
204} // End namespace Foam
205
206// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207
208#include "perfectGasI.H"
209
210#ifdef NoRepository
211 #include "perfectGas.C"
212#endif
213
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215
216#endif
217
218// ************************************************************************* //
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: perfectGas.H:91
static autoPtr< perfectGas > New(const dictionary &dict)
Definition: perfectGasI.H:65
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
Definition: perfectGasI.H:97
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: perfectGasI.H:83
static word typeName()
Return the instantiated type name.
Definition: perfectGas.H:116
autoPtr< perfectGas > clone() const
Construct and return a clone.
Definition: perfectGasI.H:57
perfectGas(const word &name, const perfectGas &)
Construct as named copy.
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
Definition: perfectGasI.H:111
void operator+=(const perfectGas &)
Definition: perfectGasI.H:141
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: perfectGasI.H:132
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: perfectGas.H:128
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: perfectGas.H:125
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectGasI.H:125
void operator*=(const scalar)
Definition: perfectGasI.H:148
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()
dictionary dict