Boussinesq.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) 2015-2021 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::Boussinesq
26 
27 Description
28  Incompressible gas equation of state using the Boussinesq approximation for
29  the density as a function of temperature only:
30 
31  \verbatim
32  rho = rho0*(1 - beta*(T - T0))
33  \endverbatim
34 
35  Coefficient mixing is very inaccurate and not supported,
36  so this equation of state is not applicable to mixtures.
37 
38 Usage
39  \table
40  Property | Description
41  rho0 | Reference density
42  T0 | Reference temperature
43  beta | Coefficient of thermal expansion
44  \endtable
45 
46  Example specification of the Boussinesq equation of state:
47  \verbatim
48  equationOfState
49  {
50  rho0 1;
51  T0 300;
52  beta 3e-03;
53  }
54  \endverbatim
55 
56 SourceFiles
57  BoussinesqI.H
58  Boussinesq.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef Boussinesq_H
63 #define Boussinesq_H
64 
65 #include "autoPtr.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of friend functions and operators
73 
74 template<class Specie> class Boussinesq;
75 
76 template<class Specie>
77 inline Boussinesq<Specie> operator+
78 (
79  const Boussinesq<Specie>&,
80  const Boussinesq<Specie>&
81 );
82 
83 template<class Specie>
84 inline Boussinesq<Specie> operator*
85 (
86  const scalar,
87  const Boussinesq<Specie>&
88 );
89 
90 template<class Specie>
91 inline Boussinesq<Specie> operator==
92 (
93  const Boussinesq<Specie>&,
94  const Boussinesq<Specie>&
95 );
96 
97 template<class Specie>
98 Ostream& operator<<
99 (
100  Ostream&,
101  const Boussinesq<Specie>&
102 );
103 
104 
105 /*---------------------------------------------------------------------------*\
106  Class Boussinesq Declaration
107 \*---------------------------------------------------------------------------*/
108 
109 template<class Specie>
110 class Boussinesq
111 :
112  public Specie
113 {
114  // Private Data
115 
116  //- Reference density
117  scalar rho0_;
118 
119  //- Reference temperature
120  scalar T0_;
121 
122  //- Thermal expansion coefficient
123  scalar beta_;
124 
125 
126 public:
127 
128  // Constructors
129 
130  //- Construct from components
131  inline Boussinesq
132  (
133  const Specie& sp,
134  const scalar rho0,
135  const scalar T0,
136  const scalar beta
137  );
138 
139  //- Construct from dictionary
140  Boussinesq(const dictionary& dict);
141 
142  //- Construct as named copy
143  inline Boussinesq
144  (
145  const word& name,
146  const Boussinesq&
147  );
148 
149  //- Construct and return a clone
150  inline autoPtr<Boussinesq> clone() const;
151 
152  // Selector from dictionary
153  inline static autoPtr<Boussinesq> New
154  (
155  const dictionary& dict
156  );
157 
158 
159  // Member Functions
160 
161  //- Return the instantiated type name
162  static word typeName()
163  {
164  return
165  "Boussinesq<"
166  + word(Specie::typeName_()) + '>';
167  }
168 
169 
170  // Fundamental properties
171 
172  //- Is the equation of state is incompressible i.e. rho != f(p)
173  static const bool incompressible = true;
174 
175  //- Is the equation of state is isochoric i.e. rho = const
176  static const bool isochoric = false;
177 
178  //- Return density [kg/m^3]
179  inline scalar rho(scalar p, scalar T) const;
180 
181  //- Return enthalpy contribution [J/kg]
182  inline scalar H(const scalar p, const scalar T) const;
183 
184  //- Return Cp contribution [J/(kg K]
185  inline scalar Cp(scalar p, scalar T) const;
186 
187  //- Return internal energy contribution [J/kg]
188  inline scalar E(const scalar p, const scalar T) const;
189 
190  //- Return Cv contribution [J/(kg K]
191  inline scalar Cv(scalar p, scalar T) const;
192 
193  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
194  inline scalar Sp(const scalar p, const scalar T) const;
195 
196  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
197  inline scalar Sv(const scalar p, const scalar T) const;
198 
199  //- Return compressibility [s^2/m^2]
200  inline scalar psi(scalar p, scalar T) const;
201 
202  //- Return compression factor []
203  inline scalar Z(scalar p, scalar T) const;
204 
205  //- Return (Cp - Cv) [J/(kg K]
206  inline scalar CpMCv(scalar p, scalar T) const;
207 
208  //- Return volumetric coefficient of thermal expansion [1/T]
209  inline scalar alphav(const scalar p, const scalar T) const;
210 
211 
212  // IO
213 
214  //- Write to Ostream
215  void write(Ostream& os) const;
216 
217 
218  // Member Operators
219 
220  inline void operator+=(const Boussinesq&);
221  inline void operator*=(const scalar);
222 
223 
224  // Friend operators
225 
226  friend Boussinesq operator+ <Specie>
227  (
228  const Boussinesq&,
229  const Boussinesq&
230  );
231 
232  friend Boussinesq operator* <Specie>
233  (
234  const scalar s,
235  const Boussinesq&
236  );
237 
238  friend Boussinesq operator== <Specie>
239  (
240  const Boussinesq&,
241  const Boussinesq&
242  );
243 
244 
245  // Ostream Operator
246 
247  friend Ostream& operator<< <Specie>
248  (
249  Ostream&,
250  const Boussinesq&
251  );
252 };
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace Foam
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #include "BoussinesqI.H"
262 
263 #ifdef NoRepository
264  #include "Boussinesq.C"
265 #endif
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #endif
270 
271 // ************************************************************************* //
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: Boussinesq.H:184
dictionary dict
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: BoussinesqI.H:137
Boussinesq(const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
Construct from components.
Definition: BoussinesqI.H:33
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: BoussinesqI.H:87
scalar rho0
static autoPtr< Boussinesq > New(const dictionary &dict)
Definition: BoussinesqI.H:72
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Definition: BoussinesqI.H:60
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: BoussinesqI.H:104
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: BoussinesqI.H:118
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 psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: BoussinesqI.H:149
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:160
A class for handling words, derived from string.
Definition: word.H:59
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: BoussinesqI.H:97
void operator+=(const Boussinesq &)
Definition: BoussinesqI.H:195
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void write(Ostream &os) const
Write to Ostream.
Definition: Boussinesq.C:47
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static word typeName()
Return the instantiated type name.
Definition: Boussinesq.H:173
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: Boussinesq.H:187
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: BoussinesqI.H:126
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: BoussinesqI.H:171
void operator*=(const scalar)
Definition: BoussinesqI.H:204
volScalarField & p
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: BoussinesqI.H:111
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:85
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: BoussinesqI.H:182
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22