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-2023 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,
35  const scalar rho0,
36  const scalar T0,
37  const scalar beta
38 )
39 :
40  Specie(sp),
41  rho0_(rho0),
42  T0_(T0),
43  beta_(beta)
44 {}
45 
46 
47 template<class Specie>
49 (
50  const word& name,
51  const Boussinesq<Specie>& b
52 )
53 :
54  Specie(name, b),
55  rho0_(b.rho0_),
56  T0_(b.T0_),
57  beta_(b.beta_)
58 {}
59 
60 
61 template<class Specie>
64 {
66  (
67  new Boussinesq<Specie>(*this)
68  );
69 }
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class Specie>
75 inline Foam::scalar Foam::Boussinesq<Specie>::rho
76 (
77  scalar p,
78  scalar T
79 ) const
80 {
81  return rho0_*(1.0 - beta_*(T - T0_));
82 }
83 
84 
85 template<class Specie>
86 inline Foam::scalar Foam::Boussinesq<Specie>::H(scalar p, scalar T) const
87 {
88  return p/this->rho(p, T);
89 }
90 
91 
92 template<class Specie>
93 inline Foam::scalar Foam::Boussinesq<Specie>::Cp(scalar p, scalar T) const
94 {
95  return 0;
96 }
97 
98 
99 template<class Specie>
100 inline Foam::scalar Foam::Boussinesq<Specie>::E(scalar p, scalar T) const
101 {
102  return 0;
103 }
104 
105 
106 template<class Specie>
107 inline Foam::scalar Foam::Boussinesq<Specie>::Cv(scalar p, scalar T) const
108 {
109  return 0;
110 }
111 
112 
113 template<class Specie>
114 inline Foam::scalar Foam::Boussinesq<Specie>::Sp
115 (
116  scalar p,
117  scalar T
118 ) const
119 {
120  return 0;
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::Boussinesq<Specie>::Sv
126 (
127  scalar p,
128  scalar T
129 ) const
130 {
132  return 0;
133 }
134 
135 
136 template<class Specie>
137 inline Foam::scalar Foam::Boussinesq<Specie>::psi
138 (
139  scalar p,
140  scalar T
141 ) const
142 {
143  return 0;
144 }
145 
146 
147 template<class Specie>
148 inline Foam::scalar Foam::Boussinesq<Specie>::Z
149 (
150  scalar p,
151  scalar T
152 ) const
153 {
154  return p/(rho(p, T)*this->R()*T);
155 }
156 
157 
158 template<class Specie>
160 (
161  scalar p,
162  scalar T
163 ) const
164 {
165  return 0;
166 }
167 
168 
169 template<class Specie>
171 (
172  scalar p,
173  scalar T
174 ) const
175 {
176  return rho0_/this->rho(p, T)*beta_;
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
181 
182 template<class Specie>
184 (
185  const Boussinesq<Specie>& b
186 )
187 {
189 }
190 
191 
192 template<class Specie>
193 inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
194 {
195  Specie::operator*=(s);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
200 
201 template<class Specie>
202 inline Foam::Boussinesq<Specie> Foam::operator+
203 (
204  const Boussinesq<Specie>& b1,
205  const Boussinesq<Specie>& b2
206 )
207 {
209  return b1;
210 }
211 
212 
213 template<class Specie>
214 inline Foam::Boussinesq<Specie> Foam::operator*
215 (
216  const scalar s,
217  const Boussinesq<Specie>& b
218 )
219 {
220  return Boussinesq<Specie>
221  (
222  s*static_cast<const Specie&>(b),
223  b.rho0_,
224  b.T0_,
225  b.beta_
226  );
227 }
228 
229 
230 template<class Specie>
231 inline Foam::Boussinesq<Specie> Foam::operator==
232 (
233  const Boussinesq<Specie>& b1,
234  const Boussinesq<Specie>& b2
235 )
236 {
237  noCoefficientMixing(Boussinesq);
238  return b1;
239 }
240 
241 
242 // ************************************************************************* //
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:124
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: BoussinesqI.H:107
autoPtr< Boussinesq > clone() const
Construct and return a clone.
Definition: BoussinesqI.H:63
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: BoussinesqI.H:126
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: BoussinesqI.H:100
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: BoussinesqI.H:138
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: BoussinesqI.H:86
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: BoussinesqI.H:171
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: BoussinesqI.H:76
Boussinesq(const Specie &sp, const scalar rho0, const scalar T0, const scalar beta)
Construct from components.
Definition: BoussinesqI.H:33
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: BoussinesqI.H:160
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: BoussinesqI.H:93
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: BoussinesqI.H:115
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: BoussinesqI.H:149
void operator*=(const scalar)
Definition: BoussinesqI.H:193
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho0
volScalarField & p
scalar T0
Definition: createFields.H:22
#define noCoefficientMixing(Type)
Definition: specie.H:159