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