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-2017 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 0;
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>::S
112 (
113  scalar p,
114  scalar T
115 ) const
116 {
117  return 0;
118 }
119 
120 
121 template<class Specie>
122 inline Foam::scalar Foam::Boussinesq<Specie>::psi
123 (
124  scalar p,
125  scalar T
126 ) const
127 {
128  return 0;
129 }
130 
131 
132 template<class Specie>
133 inline Foam::scalar Foam::Boussinesq<Specie>::Z
134 (
135  scalar p,
136  scalar T
137 ) const
138 {
139  return 0;
140 }
141 
142 
143 template<class Specie>
144 inline Foam::scalar Foam::Boussinesq<Specie>::CpMCv
145 (
146  scalar p,
147  scalar T
148 ) const
149 {
150  return this->R();
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
155 
156 template<class Specie>
157 inline void Foam::Boussinesq<Specie>::operator=
158 (
159  const Boussinesq<Specie>& b
160 )
161 {
162  Specie::operator=(b);
163 
164  rho0_ = b.rho0_;
165  T0_ = b.T0_;
166  beta_ = b.beta_;
167 }
168 
169 
170 template<class Specie>
171 inline void Foam::Boussinesq<Specie>::operator+=
172 (
173  const Boussinesq<Specie>& b
174 )
175 {
176  scalar Y1 = this->Y();
177  Specie::operator+=(b);
178 
179  if (mag(this->Y()) > SMALL)
180  {
181  Y1 /= this->Y();
182  const scalar Y2 = b.Y()/this->Y();
183 
184  rho0_ = Y1*rho0_ + Y2*b.rho0_;
185  T0_ = Y1*T0_ + Y2*b.T0_;
186  beta_ = Y1*beta_ + Y2*b.beta_;
187  }
188 }
189 
190 
191 template<class Specie>
192 inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
193 {
194  Specie::operator*=(s);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
199 
200 template<class Specie>
201 inline Foam::Boussinesq<Specie> Foam::operator+
202 (
203  const Boussinesq<Specie>& b1,
204  const Boussinesq<Specie>& b2
205 )
206 {
207  Specie sp(static_cast<const Specie&>(b1) + static_cast<const Specie&>(b2));
208 
209  if (mag(sp.Y()) < SMALL)
210  {
211  return Boussinesq<Specie>
212  (
213  sp,
214  b1.rho0_,
215  b1.T0_,
216  b1.beta_
217  );
218  }
219  else
220  {
221  const scalar Y1 = b1.Y()/sp.Y();
222  const scalar Y2 = b2.Y()/sp.Y();
223 
224  return Boussinesq<Specie>
225  (
226  sp,
227  Y1*b1.rho0_ + Y2*b2.rho0_,
228  Y1*b1.T0_ + Y2*b2.T0_,
229  Y1*b1.beta_ + Y2*b2.beta_
230  );
231  }
232 }
233 
234 
235 template<class Specie>
236 inline Foam::Boussinesq<Specie> Foam::operator*
237 (
238  const scalar s,
239  const Boussinesq<Specie>& b
240 )
241 {
242  return Boussinesq<Specie>
243  (
244  s*static_cast<const Specie&>(b),
245  b.rho0_,
246  b.T0_,
247  b.beta_
248  );
249 }
250 
251 
252 template<class Specie>
253 inline Foam::Boussinesq<Specie> Foam::operator==
254 (
255  const Boussinesq<Specie>& b1,
256  const Boussinesq<Specie>& b2
257 )
258 {
259  Specie sp(static_cast<const Specie&>(b1) == static_cast<const Specie&>(b2));
260 
261  const scalar Y1 = b1.Y()/sp.Y();
262  const scalar Y2 = b2.Y()/sp.Y();
263 
264  return Boussinesq<Specie>
265  (
266  sp,
267  Y2*b2.rho0_ - Y1*b1.rho0_,
268  Y2*b2.T0_ - Y1*b1.T0_,
269  Y2*b2.beta_ - Y1*b1.beta_
270  );
271 }
272 
273 
274 // ************************************************************************* //
PtrList< volScalarField > & Y1
Definition: YEqns.H:8
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:137
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
Definition: BoussinesqI.H:112
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
PtrList< volScalarField > & Y2
Definition: YEqns.H:9
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
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:123
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:134
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
#define R(A, B, C, D, E, F, K, M)
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:145
void operator*=(const scalar)
Definition: BoussinesqI.H:192
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:52