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-2016 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>::h(scalar p, scalar T) const
125 {
126  return 0;
127 }
128 
129 
130 template<class Specie>
131 inline Foam::scalar Foam::Boussinesq<Specie>::cp(scalar p, scalar T) const
132 {
133  return 0;
134 }
135 
136 
137 template<class Specie>
138 inline Foam::scalar Foam::Boussinesq<Specie>::s
139 (
140  scalar p,
141  scalar T
142 ) const
143 {
144  return 0;
145 }
146 
147 
148 template<class Specie>
149 inline Foam::scalar Foam::Boussinesq<Specie>::psi
150 (
151  scalar p,
152  scalar T
153 ) const
154 {
155  return 0;
156 }
157 
158 
159 template<class Specie>
160 inline Foam::scalar Foam::Boussinesq<Specie>::Z
161 (
162  scalar p,
163  scalar T
164 ) const
165 {
166  return 0;
167 }
168 
169 
170 template<class Specie>
171 inline Foam::scalar Foam::Boussinesq<Specie>::cpMcv
172 (
173  scalar p,
174  scalar T
175 ) const
176 {
177  return RR;
178 }
179 
180 
181 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
182 
183 template<class Specie>
184 inline void Foam::Boussinesq<Specie>::operator=
185 (
186  const Boussinesq<Specie>& b
187 )
188 {
189  Specie::operator=(b);
190 
191  rho0_ = b.rho0_;
192  T0_ = b.T0_;
193  beta_ = b.beta_;
194 }
195 
196 
197 template<class Specie>
198 inline void Foam::Boussinesq<Specie>::operator+=
199 (
200  const Boussinesq<Specie>& b
201 )
202 {
203  scalar molr1 = this->nMoles();
204  Specie::operator+=(b);
205  molr1 /= this->nMoles();
206  scalar molr2 = b.nMoles()/this->nMoles();
207 
208  rho0_ = molr1*rho0_ + molr2*b.rho0_;
209  T0_ = molr1*T0_ + molr2*b.T0_;
210  beta_ = molr1*beta_ + molr2*b.beta_;
211 }
212 
213 
214 template<class Specie>
215 inline void Foam::Boussinesq<Specie>::operator-=
216 (
217  const Boussinesq<Specie>& b
218 )
219 {
220  Specie::operator-=(b);
221  rho0_ = b.rho0_;
222  T0_ = b.T0_;
223  beta_ = b.beta_;
224 }
225 
226 
227 template<class Specie>
228 inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
229 {
230  Specie::operator*=(s);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
235 
236 template<class Specie>
237 inline Foam::Boussinesq<Specie> Foam::operator+
238 (
239  const Boussinesq<Specie>& b1,
240  const Boussinesq<Specie>& b2
241 )
242 {
243  scalar nMoles = b1.nMoles() + b2.nMoles();
244  scalar molr1 = b1.nMoles()/nMoles;
245  scalar molr2 = b2.nMoles()/nMoles;
246 
247  return Boussinesq<Specie>
248  (
249  static_cast<const Specie&>(b1)
250  + static_cast<const Specie&>(b2),
251  molr1*b1.rho0_ + molr2*b2.rho0_,
252  molr1*b1.T0_ + molr2*b2.T0_,
253  molr1*b1.beta_ + molr2*b2.beta_
254  );
255 }
256 
257 
258 template<class Specie>
259 inline Foam::Boussinesq<Specie> Foam::operator-
260 (
261  const Boussinesq<Specie>& b1,
262  const Boussinesq<Specie>& b2
263 )
264 {
265  return Boussinesq<Specie>
266  (
267  static_cast<const Specie&>(b1)
268  - static_cast<const Specie&>(b2),
269  b1.rho0_ - b2.rho0_,
270  b1.T0_ - b2.T0_,
271  b1.beta_ - b2.beta_
272  );
273 }
274 
275 
276 template<class Specie>
277 inline Foam::Boussinesq<Specie> Foam::operator*
278 (
279  const scalar s,
280  const Boussinesq<Specie>& b
281 )
282 {
283  return Boussinesq<Specie>
284  (
285  s*static_cast<const Specie&>(b),
286  b.rho0_,
287  b.T0_,
288  b.beta_
289  );
290 }
291 
292 
293 template<class Specie>
294 inline Foam::Boussinesq<Specie> Foam::operator==
295 (
296  const Boussinesq<Specie>& pg1,
297  const Boussinesq<Specie>& pg2
298 )
299 {
300  return pg2 - pg1;
301 }
302 
303 
304 // ************************************************************************* //
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Definition: BoussinesqI.H:73
dictionary dict
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Definition: BoussinesqI.H:139
scalar h(const scalar p, const scalar T) const
Return enthalpy departure [J/kmol].
Definition: BoussinesqI.H:124
static autoPtr< Boussinesq > New(Istream &is)
Definition: BoussinesqI.H:85
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:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:161
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))
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Definition: BoussinesqI.H:150
scalar cp(scalar p, scalar T) const
Return cp departure [J/(kmol K].
Definition: BoussinesqI.H:131
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
Definition: BoussinesqI.H:172
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: BoussinesqI.H:114
const scalar RR
Universal gas constant (default in [J/(kmol K)])
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
void operator*=(const scalar)
Definition: BoussinesqI.H:228
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:52