perfectFluidI.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) 2012-2019 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 "perfectFluid.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const scalar R,
35  const scalar rho0
36 )
37 :
38  Specie(sp),
39  R_(R),
40  rho0_(rho0)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Specie>
48 (
49  const word& name,
50  const perfectFluid<Specie>& pf
51 )
52 :
53  Specie(name, pf),
54  R_(pf.R_),
55  rho0_(pf.rho0_)
56 {}
57 
58 
59 template<class Specie>
62 {
64 }
65 
66 
67 template<class Specie>
70 (
71  const dictionary& dict
72 )
73 {
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class Specie>
81 inline Foam::scalar Foam::perfectFluid<Specie>::R() const
82 {
83  return R_;
84 }
85 
86 
87 template<class Specie>
88 inline Foam::scalar Foam::perfectFluid<Specie>::rho(scalar p, scalar T) const
89 {
90  return rho0_ + p/(this->R()*T);
91 }
92 
93 
94 template<class Specie>
95 inline Foam::scalar Foam::perfectFluid<Specie>::H(scalar p, scalar T) const
96 {
97  return p/rho(p, T) - Pstd/rho(Pstd, T);
98 }
99 
100 
101 template<class Specie>
102 inline Foam::scalar Foam::perfectFluid<Specie>::Cp(scalar p, scalar T) const
103 {
104  return
105  sqr(p/(rho(p, T)*T))/this->R()
106  - sqr(Pstd/(rho(Pstd, T)*T))/this->R();
107 }
108 
109 
110 template<class Specie>
111 inline Foam::scalar Foam::perfectFluid<Specie>::E(scalar p, scalar T) const
112 {
113  return 0;
114 }
115 
116 
117 template<class Specie>
118 inline Foam::scalar Foam::perfectFluid<Specie>::Cv(scalar p, scalar T) const
119 {
120  return 0;
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::perfectFluid<Specie>::S(scalar p, scalar T) const
126 {
127  return -this->R()*log(p/Pstd);
128 }
129 
130 
131 template<class Specie>
132 inline Foam::scalar Foam::perfectFluid<Specie>::psi(scalar p, scalar T) const
133 {
134  return 1.0/(this->R()*T);
135 }
136 
137 
138 template<class Specie>
139 inline Foam::scalar Foam::perfectFluid<Specie>::Z(scalar p, scalar T) const
140 {
141  return 1;
142 }
143 
144 
145 template<class Specie>
146 inline Foam::scalar Foam::perfectFluid<Specie>::CpMCv(scalar p, scalar T) const
147 {
148  const scalar R = this->R();
149  const scalar rho = this->rho(p, T);
150 
151  return R*sqr(p/(rho*R*T));
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
156 
157 template<class Specie>
158 inline void Foam::perfectFluid<Specie>::operator+=
159 (
160  const perfectFluid<Specie>& pf
161 )
162 {
163  scalar Y1 = this->Y();
164  Specie::operator+=(pf);
165 
166  if (mag(this->Y()) > small)
167  {
168  Y1 /= this->Y();
169  const scalar Y2 = pf.Y()/this->Y();
170 
171  R_ = 1.0/(Y1/R_ + Y2/pf.R_);
172  rho0_ = Y1*rho0_ + Y2*pf.rho0_;
173  }
174 }
175 
176 
177 template<class Specie>
178 inline void Foam::perfectFluid<Specie>::operator*=(const scalar s)
179 {
180  Specie::operator*=(s);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
185 
186 template<class Specie>
187 inline Foam::perfectFluid<Specie> Foam::operator+
188 (
189  const perfectFluid<Specie>& pf1,
190  const perfectFluid<Specie>& pf2
191 )
192 {
193  Specie sp
194  (
195  static_cast<const Specie&>(pf1)
196  + static_cast<const Specie&>(pf2)
197  );
198 
199  if (mag(sp.Y()) < small)
200  {
201  return perfectFluid<Specie>
202  (
203  sp,
204  pf1.R_,
205  pf1.rho0_
206  );
207  }
208  else
209  {
210  const scalar Y1 = pf1.Y()/sp.Y();
211  const scalar Y2 = pf2.Y()/sp.Y();
212 
213  return perfectFluid<Specie>
214  (
215  sp,
216  1.0/(Y1/pf1.R_ + Y2/pf2.R_),
217  Y1*pf1.rho0_ + Y2*pf2.rho0_
218  );
219  }
220 }
221 
222 
223 template<class Specie>
224 inline Foam::perfectFluid<Specie> Foam::operator*
225 (
226  const scalar s,
227  const perfectFluid<Specie>& pf
228 )
229 {
230  return perfectFluid<Specie>
231  (
232  s*static_cast<const Specie&>(pf),
233  pf.R_,
234  pf.rho0_
235  );
236 }
237 
238 
239 template<class Specie>
240 inline Foam::perfectFluid<Specie> Foam::operator==
241 (
242  const perfectFluid<Specie>& pf1,
243  const perfectFluid<Specie>& pf2
244 )
245 {
246  Specie sp
247  (
248  static_cast<const Specie&>(pf1)
249  == static_cast<const Specie&>(pf2)
250  );
251 
252  if (mag(sp.Y()) < small)
253  {
254  return perfectFluid<Specie>
255  (
256  sp,
257  pf1.R_,
258  pf1.rho0_
259  );
260  }
261  else
262  {
263  const scalar Y1 = pf1.Y()/sp.Y();
264  const scalar Y2 = pf2.Y()/sp.Y();
265  const scalar oneByR = Y2/pf2.R_ - Y1/pf1.R_;
266 
267  return perfectFluid<Specie>
268  (
269  sp,
270  mag(oneByR) < small ? great : 1/oneByR,
271  Y2*pf2.rho0_ - Y1*pf1.rho0_
272  );
273  }
274 }
275 
276 
277 // ************************************************************************* //
autoPtr< perfectFluid > clone() const
Construct and return a clone.
Definition: perfectFluidI.H:61
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensionedSymmTensor sqr(const dimensionedVector &dv)
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
Construct from components.
Definition: perfectFluidI.H:32
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectFluidI.H:88
static autoPtr< perfectFluid > New(const dictionary &dict)
Definition: perfectFluidI.H:70
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar R() const
Return fluid constant [J/kg/K].
Definition: perfectFluidI.H:81
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))
A class for handling words, derived from string.
Definition: word.H:59
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
void operator*=(const scalar)
const dimensionedScalar Pstd
Standard pressure.
scalar Z(scalar p, scalar T) const
Return compression factor [].
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: perfectFluidI.H:95
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
const volScalarField & T
Perfect gas equation of state.
Definition: perfectFluid.H:47
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
#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 Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
volScalarField & p
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].