Boussinesq.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) 2015-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::Boussinesq
28
29Group
30 grpSpecieEquationOfState
31
32Description
33 Incompressible gas equation of state using the Boussinesq approximation for
34 the density as a function of temperature only:
35
36 \verbatim
37 rho = rho0*(1 - beta*(T - T0))
38 \endverbatim
39
40SourceFiles
41 BoussinesqI.H
42 Boussinesq.C
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef Boussinesq_H
47#define Boussinesq_H
48
49#include "autoPtr.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56// Forward declaration of friend functions and operators
57
58template<class Specie> class Boussinesq;
59
60template<class Specie>
61inline Boussinesq<Specie> operator+
62(
63 const Boussinesq<Specie>&,
65);
66
67template<class Specie>
68inline Boussinesq<Specie> operator*
69(
70 const scalar,
72);
73
74template<class Specie>
75inline Boussinesq<Specie> operator==
76(
77 const Boussinesq<Specie>&,
79);
80
81template<class Specie>
82Ostream& operator<<
83(
84 Ostream&,
86);
87
88
89/*---------------------------------------------------------------------------*\
90 Class Boussinesq Declaration
91\*---------------------------------------------------------------------------*/
92
93template<class Specie>
94class Boussinesq
95:
96 public Specie
97{
98 // Private data
99
100 //- Reference density
101 scalar rho0_;
102
103 //- Reference temperature
104 scalar T0_;
105
106 //- Thermal expansion coefficient
107 scalar beta_;
108
109
110public:
111
112 // Constructors
113
114 //- Construct from components
115 inline Boussinesq
116 (
117 const Specie& sp,
118 const scalar rho0,
119 const scalar T0,
120 const scalar beta
121 );
122
123 //- Construct from dictionary
124 Boussinesq(const dictionary& dict);
125
126 //- Construct as named copy
127 inline Boussinesq
128 (
129 const word& name,
130 const Boussinesq&
131 );
132
133 //- Construct and return a clone
134 inline autoPtr<Boussinesq> clone() const;
135
136 // Selector from dictionary
137 inline static autoPtr<Boussinesq> New
138 (
139 const dictionary& dict
140 );
141
142
143 // Member functions
144
145 //- Return the instantiated type name
146 static word typeName()
147 {
148 return
149 "Boussinesq<"
150 + word(Specie::typeName_()) + '>';
151 }
152
153
154 // Fundamental properties
155
156 //- Is the equation of state is incompressible i.e. rho != f(p)
157 static const bool incompressible = true;
158
159 //- Is the equation of state is isochoric i.e. rho = const
160 static const bool isochoric = false;
161
162 //- Return density [kg/m^3]
163 inline scalar rho(scalar p, scalar T) const;
164
165 //- Return enthalpy departure [J/kg]
166 inline scalar H(const scalar p, const scalar T) const;
167
168 //- Return Cp departure [J/(kg K]
169 inline scalar Cp(scalar p, scalar T) const;
170
171 //- Return internal energy departure [J/kg]
172 inline scalar E(const scalar p, const scalar T) const;
173
174 //- Return Cv departure [J/(kg K]
175 inline scalar Cv(scalar p, scalar T) const;
176
177 //- Return entropy [J/(kg K)]
178 inline scalar S(const scalar p, const scalar T) const;
179
180 //- Return compressibility rho/p [s^2/m^2]
181 inline scalar psi(scalar p, scalar T) const;
182
183 //- Return compression factor []
184 inline scalar Z(scalar p, scalar T) const;
185
186 //- Return (Cp - Cv) [J/(kg K]
187 inline scalar CpMCv(scalar p, scalar T) const;
188
189
190 // IO
191
192 //- Write to Ostream
193 void write(Ostream& os) const;
194
195
196 // Member operators
197
198 inline void operator=(const Boussinesq&);
199 inline void operator+=(const Boussinesq&);
200 inline void operator*=(const scalar);
201
202
203 // Friend operators
205 friend Boussinesq operator+ <Specie>
206 (
207 const Boussinesq&,
208 const Boussinesq&
209 );
211 friend Boussinesq operator* <Specie>
212 (
213 const scalar s,
214 const Boussinesq&
215 );
217 friend Boussinesq operator== <Specie>
218 (
219 const Boussinesq&,
220 const Boussinesq&
221 );
222
223
224 // Ostream Operator
226 friend Ostream& operator<< <Specie>
227 (
228 Ostream&,
229 const Boussinesq&
230 );
231};
232
233
234// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235
236} // End namespace Foam
237
238// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239
240#include "BoussinesqI.H"
241
242#ifdef NoRepository
243 #include "Boussinesq.C"
244#endif
245
246// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247
248#endif
249
250// ************************************************************************* //
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:96
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Definition: BoussinesqI.H:62
void operator+=(const Boussinesq &)
Definition: BoussinesqI.H:181
void operator=(const Boussinesq &)
Definition: BoussinesqI.H:167
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
Definition: BoussinesqI.H:106
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: BoussinesqI.H:93
static autoPtr< Boussinesq > New(const dictionary &dict)
Definition: BoussinesqI.H:71
static word typeName()
Return the instantiated type name.
Definition: Boussinesq.H:145
Boussinesq(const word &name, const Boussinesq &)
Construct as named copy.
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
Definition: BoussinesqI.H:121
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: BoussinesqI.H:154
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: Boussinesq.H:159
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: Boussinesq.H:156
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:143
void operator*=(const scalar)
Definition: BoussinesqI.H:201
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
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 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
scalar T0
Definition: createFields.H:22
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)