perfectFluidI.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) 2012-2015 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 {
73 }
74 
75 
76 template<class Specie>
79 (
80  const dictionary& dict
81 )
82 {
84 }
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
89 template<class Specie>
90 inline Foam::scalar Foam::perfectFluid<Specie>::R() const
91 {
92  return R_;
93 }
94 
95 
96 template<class Specie>
97 inline Foam::scalar Foam::perfectFluid<Specie>::rho(scalar p, scalar T) const
98 {
99  return rho0_ + p/(this->R()*T);
100 }
101 
102 
103 template<class Specie>
104 inline Foam::scalar Foam::perfectFluid<Specie>::s(scalar p, scalar T) const
105 {
106  return -RR*log(p/Pstd);
107 }
108 
109 
110 template<class Specie>
111 inline Foam::scalar Foam::perfectFluid<Specie>::psi(scalar p, scalar T) const
112 {
113  return 1.0/(this->R()*T);
114 }
115 
116 
117 template<class Specie>
118 inline Foam::scalar Foam::perfectFluid<Specie>::Z(scalar p, scalar T) const
119 {
120  return 1;
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::perfectFluid<Specie>::cpMcv(scalar p, scalar T) const
126 {
127  return 0;
128 }
129 
130 
131 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
132 
133 template<class Specie>
135 (
136  const perfectFluid<Specie>& pf
137 )
138 {
139  scalar molr1 = this->nMoles();
140 
141  Specie::operator+=(pf);
142 
143  molr1 /= this->nMoles();
144  scalar molr2 = pf.nMoles()/this->nMoles();
145 
146  R_ = 1.0/(molr1/R_ + molr2/pf.R_);
147  rho0_ = molr1*rho0_ + molr2*pf.rho0_;
148 }
149 
150 
151 template<class Specie>
153 (
154  const perfectFluid<Specie>& pf
155 )
156 {
157  scalar molr1 = this->nMoles();
158 
159  Specie::operator-=(pf);
160 
161  molr1 /= this->nMoles();
162  scalar molr2 = pf.nMoles()/this->nMoles();
163 
164  R_ = 1.0/(molr1/R_ - molr2/pf.R_);
165  rho0_ = molr1*rho0_ - molr2*pf.rho0_;
166 }
167 
168 
169 template<class Specie>
170 inline void Foam::perfectFluid<Specie>::operator*=(const scalar s)
171 {
172  Specie::operator*=(s);
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
177 
178 template<class Specie>
179 inline Foam::perfectFluid<Specie> Foam::operator+
180 (
181  const perfectFluid<Specie>& pf1,
182  const perfectFluid<Specie>& pf2
183 )
184 {
185  scalar nMoles = pf1.nMoles() + pf2.nMoles();
186  scalar molr1 = pf1.nMoles()/nMoles;
187  scalar molr2 = pf2.nMoles()/nMoles;
188 
189  return perfectFluid<Specie>
190  (
191  static_cast<const Specie&>(pf1)
192  + static_cast<const Specie&>(pf2),
193  1.0/(molr1/pf1.R_ + molr2/pf2.R_),
194  molr1*pf1.rho0_ + molr2*pf2.rho0_
195  );
196 }
197 
198 
199 template<class Specie>
200 inline Foam::perfectFluid<Specie> Foam::operator-
201 (
202  const perfectFluid<Specie>& pf1,
203  const perfectFluid<Specie>& pf2
204 )
205 {
206  scalar nMoles = pf1.nMoles() + pf2.nMoles();
207  scalar molr1 = pf1.nMoles()/nMoles;
208  scalar molr2 = pf2.nMoles()/nMoles;
209 
210  return perfectFluid<Specie>
211  (
212  static_cast<const Specie&>(pf1)
213  - static_cast<const Specie&>(pf2),
214  1.0/(molr1/pf1.R_ - molr2/pf2.R_),
215  molr1*pf1.rho0_ - molr2*pf2.rho0_
216  );
217 }
218 
219 
220 template<class Specie>
221 inline Foam::perfectFluid<Specie> Foam::operator*
222 (
223  const scalar s,
224  const perfectFluid<Specie>& pf
225 )
226 {
227  return perfectFluid<Specie>
228  (
229  s*static_cast<const Specie&>(pf),
230  pf.R_,
231  pf.rho0_
232  );
233 }
234 
235 
236 template<class Specie>
237 inline Foam::perfectFluid<Specie> Foam::operator==
238 (
239  const perfectFluid<Specie>& pf1,
240  const perfectFluid<Specie>& pf2
241 )
242 {
243  return pf2 - pf1;
244 }
245 
246 
247 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Definition: thermodynamicConstants.C:44
Foam::perfectFluid
Perfect gas equation of state.
Definition: perfectFluid.H:47
Foam::perfectFluid::cpMcv
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
Definition: perfectFluidI.H:125
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::perfectFluid::perfectFluid
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
Construct from components.
Definition: perfectFluidI.H:33
Foam::constant::standard::Pstd
const dimensionedScalar Pstd
Standard pressure.
Definition: thermodynamicConstants.C:46
Foam::perfectFluid::rho0_
scalar rho0_
The reference density.
Definition: perfectFluid.H:100
perfectFluid.H
Foam::perfectFluid::R_
scalar R_
Fluid constant.
Definition: perfectFluid.H:97
Foam::perfectFluid::New
static autoPtr< perfectFluid > New(Istream &is)
Definition: perfectFluidI.H:70
specie.H
Foam::perfectFluid::clone
autoPtr< perfectFluid > clone() const
Construct and return a clone.
Definition: perfectFluidI.H:62
Foam::perfectFluid::operator*=
void operator*=(const scalar)
Definition: perfectFluidI.H:170
R
#define R(A, B, C, D, E, F, K, M)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::perfectFluid::psi
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
Definition: perfectFluidI.H:111
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::perfectFluid::rho
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectFluidI.H:97
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
s
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))
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
rho0
scalar rho0
Definition: readInitialConditions.H:97
Foam::perfectFluid::R
scalar R() const
Return fluid constant [J/(kg K)].
Definition: perfectFluidI.H:90
Foam::perfectFluid::Z
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectFluidI.H:118
Foam::perfectFluid::s
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
Definition: perfectFluidI.H:104
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47