BoussinesqI.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 | Copyright (C) 2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "Boussinesq.H"
27 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Specie>
33 (
34  const Specie& sp, const scalar rho0, const scalar T0, const scalar beta
35 )
36 :
37  Specie(sp),
38  rho0_(rho0),
39  T0_(T0),
40  beta_(beta)
41 {}
42 
43 
44 template<class Specie>
46 (
47  const Boussinesq& b
48 )
49 :
50  Specie(b),
51  rho0_(b.rho0_),
52  T0_(b.T0_),
53  beta_(b.beta_)
54 {}
55 
56 
57 template<class Specie>
59 (
60  const word& name,
61  const Boussinesq<Specie>& b
62 )
63 :
64  Specie(name, b),
65  rho0_(b.rho0_),
66  T0_(b.T0_),
67  beta_(b.beta_)
68 {}
69 
70 
71 template<class Specie>
74 {
76  (
77  new Boussinesq<Specie>(*this)
78  );
79 }
80 
81 
82 template<class Specie>
85 (
86  Istream& is
87 )
88 {
90  (
91  new Boussinesq<Specie>(is)
92  );
93 }
94 
95 
96 template<class Specie>
99 (
100  const dictionary& dict
101 )
102 {
104  (
106  );
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 template<class Specie>
113 inline Foam::scalar Foam::Boussinesq<Specie>::rho
114 (
115  scalar p,
116  scalar T
117 ) const
118 {
119  return rho0_*(1.0 - beta_*(T - T0_));
120 }
121 
122 
123 template<class Specie>
124 inline Foam::scalar Foam::Boussinesq<Specie>::s
125 (
126  scalar p,
127  scalar T
128 ) const
129 {
130  return 0;
131 }
132 
133 
134 template<class Specie>
135 inline Foam::scalar Foam::Boussinesq<Specie>::psi
136 (
137  scalar p,
138  scalar T
139 ) const
140 {
141  return 0;
142 }
143 
144 
145 template<class Specie>
146 inline Foam::scalar Foam::Boussinesq<Specie>::Z
147 (
148  scalar p,
149  scalar T
150 ) const
151 {
152  return 0;
153 }
154 
155 
156 template<class Specie>
157 inline Foam::scalar Foam::Boussinesq<Specie>::cpMcv
158 (
159  scalar p,
160  scalar T
161 ) const
162 {
163  return RR;
164 }
165 
166 
167 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
168 
169 template<class Specie>
172 (
173  const Boussinesq<Specie>& b
174 )
175 {
176  Specie::operator=(b);
177 
178  rho0_ = b.rho0_;
179  T0_ = b.T0_;
180  beta_ = b.beta_;
181 
182  return *this;
183 }
184 
185 template<class Specie>
187 (
188  const Boussinesq<Specie>& b
189 )
190 {
191  scalar molr1 = this->nMoles();
192  Specie::operator+=(b);
193  molr1 /= this->nMoles();
194  scalar molr2 = b.nMoles()/this->nMoles();
195 
196  rho0_ = molr1*rho0_ + molr2*b.rho0_;
197  T0_ = molr1*T0_ + molr2*b.T0_;
198  beta_ = molr1*beta_ + molr2*b.beta_;
199 }
200 
201 
202 template<class Specie>
204 (
205  const Boussinesq<Specie>& b
206 )
207 {
208  Specie::operator-=(b);
209  rho0_ = b.rho0_;
210  T0_ = b.T0_;
211  beta_ = b.beta_;
212 }
213 
214 
215 template<class Specie>
216 inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
217 {
218  Specie::operator*=(s);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
223 
224 template<class Specie>
225 inline Foam::Boussinesq<Specie> Foam::operator+
226 (
227  const Boussinesq<Specie>& b1,
228  const Boussinesq<Specie>& b2
229 )
230 {
231  scalar nMoles = b1.nMoles() + b2.nMoles();
232  scalar molr1 = b1.nMoles()/nMoles;
233  scalar molr2 = b2.nMoles()/nMoles;
234 
235  return Boussinesq<Specie>
236  (
237  static_cast<const Specie&>(b1)
238  + static_cast<const Specie&>(b2),
239  molr1*b1.rho0_ + molr2*b2.rho0_,
240  molr1*b1.T0_ + molr2*b2.T0_,
241  molr1*b1.beta_ + molr2*b2.beta_
242  );
243 }
244 
245 
246 template<class Specie>
247 inline Foam::Boussinesq<Specie> Foam::operator-
248 (
249  const Boussinesq<Specie>& b1,
250  const Boussinesq<Specie>& b2
251 )
252 {
253  return Boussinesq<Specie>
254  (
255  static_cast<const Specie&>(b1)
256  - static_cast<const Specie&>(b2),
257  b1.rho0_ - b2.rho0_,
258  b1.T0_ - b2.T0_,
259  b1.beta_ - b2.beta_
260  );
261 }
262 
263 
264 template<class Specie>
265 inline Foam::Boussinesq<Specie> Foam::operator*
266 (
267  const scalar s,
268  const Boussinesq<Specie>& b
269 )
270 {
271  return Boussinesq<Specie>
272  (
273  s*static_cast<const Specie&>(b),
274  b.rho0_,
275  b.T0_,
276  b.beta_
277  );
278 }
279 
280 
281 template<class Specie>
282 inline Foam::Boussinesq<Specie> Foam::operator==
283 (
284  const Boussinesq<Specie>& pg1,
285  const Boussinesq<Specie>& pg2
286 )
287 {
288  return pg2 - pg1;
289 }
290 
291 
292 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Definition: thermodynamicConstants.C:44
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::Boussinesq::s
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Definition: BoussinesqI.H:125
p
p
Definition: pEqn.H:62
Foam::Boussinesq::clone
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Definition: BoussinesqI.H:73
Foam::Boussinesq::Boussinesq
Boussinesq(const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
Construct from components.
Definition: BoussinesqI.H:33
Foam::Boussinesq::T0_
scalar T0_
Reference temperature.
Definition: Boussinesq.H:105
Foam::Boussinesq::psi
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Definition: BoussinesqI.H:136
Foam::Boussinesq::cpMcv
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
Definition: BoussinesqI.H:158
Foam::Boussinesq::Z
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:147
Foam::Boussinesq::operator*=
void operator*=(const scalar)
Definition: BoussinesqI.H:216
Foam::Boussinesq
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:52
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::Boussinesq::rho
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: BoussinesqI.H:114
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
s
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))
Foam::Boussinesq::beta_
scalar beta_
Thermal expansion coefficient.
Definition: Boussinesq.H:108
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
rho0
scalar rho0
Definition: readInitialConditions.H:97
Foam::Boussinesq::rho0_
scalar rho0_
Reference density.
Definition: Boussinesq.H:102
Foam::Boussinesq::New
static autoPtr< Boussinesq > New(Istream &is)
Definition: BoussinesqI.H:85
Boussinesq.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47