All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BoussinesqI.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-2018 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 word& name,
48  const Boussinesq<Specie>& b
49 )
50 :
51  Specie(name, b),
52  rho0_(b.rho0_),
53  T0_(b.T0_),
54  beta_(b.beta_)
55 {}
56 
57 
58 template<class Specie>
61 {
63  (
64  new Boussinesq<Specie>(*this)
65  );
66 }
67 
68 
69 template<class Specie>
72 (
73  const dictionary& dict
74 )
75 {
77  (
79  );
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class Specie>
86 inline Foam::scalar Foam::Boussinesq<Specie>::rho
87 (
88  scalar p,
89  scalar T
90 ) const
91 {
92  return rho0_*(1.0 - beta_*(T - T0_));
93 }
94 
95 
96 template<class Specie>
97 inline Foam::scalar Foam::Boussinesq<Specie>::H(scalar p, scalar T) const
98 {
99  return p/this->rho(p, T);
100 }
101 
102 
103 template<class Specie>
104 inline Foam::scalar Foam::Boussinesq<Specie>::Cp(scalar p, scalar T) const
105 {
106  return 0;
107 }
108 
109 
110 template<class Specie>
111 inline Foam::scalar Foam::Boussinesq<Specie>::E(scalar p, scalar T) const
112 {
113  return 0;
114 }
115 
116 
117 template<class Specie>
118 inline Foam::scalar Foam::Boussinesq<Specie>::Cv(scalar p, scalar T) const
119 {
120  return 0;
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::Boussinesq<Specie>::S
126 (
127  scalar p,
128  scalar T
129 ) const
130 {
131  return 0;
132 }
133 
134 
135 template<class Specie>
136 inline Foam::scalar Foam::Boussinesq<Specie>::psi
137 (
138  scalar p,
139  scalar T
140 ) const
141 {
142  return 0;
143 }
144 
145 
146 template<class Specie>
147 inline Foam::scalar Foam::Boussinesq<Specie>::Z
148 (
149  scalar p,
150  scalar T
151 ) const
152 {
153  return 0;
154 }
155 
156 
157 template<class Specie>
158 inline Foam::scalar Foam::Boussinesq<Specie>::CpMCv
159 (
160  scalar p,
161  scalar T
162 ) const
163 {
164  return 0;
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
169 
170 template<class Specie>
171 inline void Foam::Boussinesq<Specie>::operator=
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 
183 
184 template<class Specie>
185 inline void Foam::Boussinesq<Specie>::operator+=
186 (
187  const Boussinesq<Specie>& b
188 )
189 {
190  scalar Y1 = this->Y();
191  Specie::operator+=(b);
192 
193  if (mag(this->Y()) > small)
194  {
195  Y1 /= this->Y();
196  const scalar Y2 = b.Y()/this->Y();
197 
198  rho0_ = Y1*rho0_ + Y2*b.rho0_;
199  T0_ = Y1*T0_ + Y2*b.T0_;
200  beta_ = Y1*beta_ + Y2*b.beta_;
201  }
202 }
203 
204 
205 template<class Specie>
206 inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
207 {
208  Specie::operator*=(s);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
213 
214 template<class Specie>
215 inline Foam::Boussinesq<Specie> Foam::operator+
216 (
217  const Boussinesq<Specie>& b1,
218  const Boussinesq<Specie>& b2
219 )
220 {
221  Specie sp(static_cast<const Specie&>(b1) + static_cast<const Specie&>(b2));
222 
223  if (mag(sp.Y()) < small)
224  {
225  return Boussinesq<Specie>
226  (
227  sp,
228  b1.rho0_,
229  b1.T0_,
230  b1.beta_
231  );
232  }
233  else
234  {
235  const scalar Y1 = b1.Y()/sp.Y();
236  const scalar Y2 = b2.Y()/sp.Y();
237 
238  return Boussinesq<Specie>
239  (
240  sp,
241  Y1*b1.rho0_ + Y2*b2.rho0_,
242  Y1*b1.T0_ + Y2*b2.T0_,
243  Y1*b1.beta_ + Y2*b2.beta_
244  );
245  }
246 }
247 
248 
249 template<class Specie>
250 inline Foam::Boussinesq<Specie> Foam::operator*
251 (
252  const scalar s,
253  const Boussinesq<Specie>& b
254 )
255 {
256  return Boussinesq<Specie>
257  (
258  s*static_cast<const Specie&>(b),
259  b.rho0_,
260  b.T0_,
261  b.beta_
262  );
263 }
264 
265 
266 template<class Specie>
267 inline Foam::Boussinesq<Specie> Foam::operator==
268 (
269  const Boussinesq<Specie>& b1,
270  const Boussinesq<Specie>& b2
271 )
272 {
273  Specie sp(static_cast<const Specie&>(b1) == static_cast<const Specie&>(b2));
274 
275  const scalar Y1 = b1.Y()/sp.Y();
276  const scalar Y2 = b2.Y()/sp.Y();
277 
278  return Boussinesq<Specie>
279  (
280  sp,
281  Y2*b2.rho0_ - Y1*b1.rho0_,
282  Y2*b2.T0_ - Y1*b1.T0_,
283  Y2*b2.beta_ - Y1*b1.beta_
284  );
285 }
286 
287 
288 // ************************************************************************* //
dictionary dict
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:158
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
Definition: BoussinesqI.H:126
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: BoussinesqI.H:87
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 departure [J/(kg K].
Definition: BoussinesqI.H:104
scalar Cv(scalar p, scalar T) const
Return Cv departure [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 rho/p [s^2/m^2].
Definition: BoussinesqI.H:137
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:148
A class for handling words, derived from string.
Definition: word.H:59
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: BoussinesqI.H:97
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
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:159
void operator*=(const scalar)
Definition: BoussinesqI.H:206
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
Definition: BoussinesqI.H:111
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:52