adiabaticPerfectFluidI.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) 2013-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 "adiabaticPerfectFluid.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const scalar p0,
35  const scalar rho0,
36  const scalar gamma,
37  const scalar B
38 )
39 :
40  Specie(sp),
41  p0_(p0),
42  rho0_(rho0),
43  gamma_(gamma),
44  B_(B)
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class Specie>
52 (
53  const word& name,
55 )
56 :
57  Specie(name, pf),
58  p0_(pf.p0_),
59  rho0_(pf.rho0_),
60  gamma_(pf.gamma_),
61  B_(pf.B_)
62 {}
63 
64 
65 template<class Specie>
68 {
70  (
72  );
73 }
74 
75 
76 template<class Specie>
79 {
81  (
83  );
84 }
85 
86 
87 template<class Specie>
90 (
91  const dictionary& dict
92 )
93 {
95  (
97  );
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class Specie>
105 (
106  scalar p,
107  scalar T
108 ) const
109 {
110  return rho0_*pow((p + B_)/(p0_ + B_), 1.0/gamma_);
111 }
112 
113 
114 template<class Specie>
115 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::s
116 (
117  scalar p,
118  scalar T
119 ) const
120 {
121  scalar n = 1 - 1.0/gamma_;
122  return
123  -pow(p0_ + B_, 1.0/gamma_)*(pow((p + B_), n) - pow((Pstd + B_), n))
124  /(rho0_*T*n);
125 }
126 
127 
128 template<class Specie>
130 (
131  scalar p,
132  scalar T
133 ) const
134 {
135  return
136  (rho0_/(gamma_*(p0_ + B_)))
137  *pow((p + B_)/(p0_ + B_), 1.0/gamma_ - 1.0);
138 }
139 
140 
141 template<class Specie>
142 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Z(scalar, scalar) const
143 {
144  return 1;
145 }
146 
147 
148 template<class Specie>
150 (
151  scalar p,
152  scalar T
153 ) const
154 {
155  return 0;
156 }
157 
158 
159 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
160 
161 template<class Specie>
162 inline void Foam::adiabaticPerfectFluid<Specie>::operator+=
163 (
165 )
166 {
167  scalar molr1 = this->nMoles();
168 
169  Specie::operator+=(pf);
170 
171  molr1 /= this->nMoles();
172  scalar molr2 = pf.nMoles()/this->nMoles();
173 
174  p0_ = molr1*p0_ + molr2*pf.p0_;
175  rho0_ = molr1*rho0_ + molr2*pf.rho0_;
176  gamma_ = molr1*gamma_ + molr2*pf.gamma_;
177  B_ = molr1*B_ + molr2*pf.B_;
178 }
179 
180 
181 template<class Specie>
182 inline void Foam::adiabaticPerfectFluid<Specie>::operator-=
183 (
185 )
186 {
187  scalar molr1 = this->nMoles();
188 
189  Specie::operator-=(pf);
190 
191  molr1 /= this->nMoles();
192  scalar molr2 = pf.nMoles()/this->nMoles();
193 
194  p0_ = molr1*p0_ - molr2*pf.p0_;
195  rho0_ = molr1*rho0_ - molr2*pf.rho0_;
196  gamma_ = molr1*gamma_ - molr2*pf.gamma_;
197  B_ = molr1*B_ - molr2*pf.B_;
198 }
199 
200 
201 template<class Specie>
203 {
204  Specie::operator*=(s);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
209 
210 template<class Specie>
211 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator+
212 (
215 )
216 {
217  scalar nMoles = pf1.nMoles() + pf2.nMoles();
218  scalar molr1 = pf1.nMoles()/nMoles;
219  scalar molr2 = pf2.nMoles()/nMoles;
220 
221  return rhoConst<Specie>
222  (
223  static_cast<const Specie&>(pf1)
224  + static_cast<const Specie&>(pf2),
225  molr1*pf1.p0_ + molr2*pf2.p0_,
226  molr1*pf1.rho0_ + molr2*pf2.rho0_,
227  molr1*pf1.gamma_ + molr2*pf2.gamma_,
228  molr1*pf1.B_ + molr2*pf2.B_
229  );
230 }
231 
232 
233 template<class Specie>
234 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator-
235 (
238 )
239 {
240  scalar nMoles = pf1.nMoles() + pf2.nMoles();
241  scalar molr1 = pf1.nMoles()/nMoles;
242  scalar molr2 = pf2.nMoles()/nMoles;
243 
244  return rhoConst<Specie>
245  (
246  static_cast<const Specie&>(pf1)
247  - static_cast<const Specie&>(pf2),
248  molr1*pf1.p0_ - molr2*pf2.p0_,
249  molr1*pf1.rho0_ - molr2*pf2.rho0_,
250  molr1*pf1.gamma_ - molr2*pf2.gamma_,
251  molr1*pf1.B_ - molr2*pf2.B_
252  );
253 }
254 
255 
256 template<class Specie>
257 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator*
258 (
259  const scalar s,
261 )
262 {
264  (
265  s*static_cast<const Specie&>(pf),
266  pf.p0_,
267  pf.rho0_,
268  pf.gamma_,
269  pf.B_
270  );
271 }
272 
273 
274 template<class Specie>
275 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator==
276 (
279 )
280 {
281  return pf2 - pf1;
282 }
283 
284 
285 // ************************************************************************* //
RhoConst (rho = const) of state.
Definition: rhoConst.H:47
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 ))
static autoPtr< adiabaticPerfectFluid > New(Istream &is)
AdiabaticPerfect gas equation of state.
scalar Z(scalar p, scalar T) const
Return compression factor [].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
A class for handling words, derived from string.
Definition: word.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dimensionedScalar Pstd
Standard pressure.
label n
dictionary dict
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar s(const scalar p, const scalar T) const
Return entropy [J/(kmol K)].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
scalar cpMcv(scalar p, scalar T) const
Return (cp - cv) [J/(kmol K].
adiabaticPerfectFluid(const Specie &sp, const scalar p0, const scalar rho0, const scalar gamma, const scalar B)
Construct from components.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117