adiabaticPerfectFluidI.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) 2013-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 "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 (
80  const dictionary& dict
81 )
82 {
84  (
86  );
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
92 template<class Specie>
94 (
95  scalar p,
96  scalar T
97 ) const
98 {
99  return rho0_*pow((p + B_)/(p0_ + B_), 1.0/gamma_);
100 }
101 
102 
103 template<class Specie>
104 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::H
105 (
106  scalar p,
107  scalar T
108 ) const
109 {
110  return 0;
111 }
112 
113 
114 template<class Specie>
115 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Cp
116 (
117  scalar p,
118  scalar T
119 ) const
120 {
121  return 0;
122 }
123 
124 
125 template<class Specie>
126 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::S
127 (
128  scalar p,
129  scalar T
130 ) const
131 {
132  scalar n = 1 - 1.0/gamma_;
133  return
134  -pow(p0_ + B_, 1.0/gamma_)*(pow((p + B_), n) - pow((Pstd + B_), n))
135  /(rho0_*T*n);
136 }
137 
138 
139 template<class Specie>
141 (
142  scalar p,
143  scalar T
144 ) const
145 {
146  return
147  (rho0_/(gamma_*(p0_ + B_)))
148  *pow((p + B_)/(p0_ + B_), 1.0/gamma_ - 1.0);
149 }
150 
151 
152 template<class Specie>
153 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Z(scalar, scalar) const
154 {
155  return 1;
156 }
157 
158 
159 template<class Specie>
161 (
162  scalar p,
163  scalar T
164 ) const
165 {
166  return 0;
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
171 
172 template<class Specie>
173 inline void Foam::adiabaticPerfectFluid<Specie>::operator+=
174 (
176 )
177 {
178  scalar Y1 = this->Y();
179  Specie::operator+=(pf);
180 
181  if (mag(this->Y()) > small)
182  {
183  Y1 /= this->Y();
184  const scalar Y2 = pf.Y()/this->Y();
185 
186  p0_ = Y1*p0_ + Y2*pf.p0_;
187  rho0_ = Y1*rho0_ + Y2*pf.rho0_;
188  gamma_ = Y1*gamma_ + Y2*pf.gamma_;
189  B_ = Y1*B_ + Y2*pf.B_;
190  }
191 }
192 
193 
194 template<class Specie>
196 {
197  Specie::operator*=(s);
198 }
199 
200 
201 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
202 
203 template<class Specie>
204 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator+
205 (
208 )
209 {
210  Specie sp
211  (
212  static_cast<const Specie&>(pf1)
213  + static_cast<const Specie&>(pf2)
214  );
215 
216  if (mag(sp.Y()) < small)
217  {
219  (
220  sp,
221  pf1.p0_,
222  pf1.rho0_,
223  pf1.gamma_,
224  pf1.B_
225  );
226  }
227  else
228  {
229  const scalar Y1 = pf1.Y()/sp.Y();
230  const scalar Y2 = pf2.Y()/sp.Y();
231 
233  (
234  sp,
235  Y1*pf1.p0_ + Y2*pf2.p0_,
236  Y1*pf1.rho0_ + Y2*pf2.rho0_,
237  Y1*pf1.gamma_ + Y2*pf2.gamma_,
238  Y1*pf1.B_ + Y2*pf2.B_
239  );
240  }
241 }
242 
243 
244 template<class Specie>
245 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator*
246 (
247  const scalar s,
249 )
250 {
252  (
253  s*static_cast<const Specie&>(pf),
254  pf.p0_,
255  pf.rho0_,
256  pf.gamma_,
257  pf.B_
258  );
259 }
260 
261 
262 template<class Specie>
263 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator==
264 (
267 )
268 {
269  Specie sp
270  (
271  static_cast<const Specie&>(pf1)
272  == static_cast<const Specie&>(pf2)
273  );
274 
275  const scalar Y1 = pf1.Y()/sp.Y();
276  const scalar Y2 = pf2.Y()/sp.Y();
277 
279  (
280  sp,
281  Y2*pf2.p0_ - Y1*pf1.p0_,
282  Y2*pf2.rho0_ - Y1*pf1.rho0_,
283  Y2*pf2.gamma_ - Y1*pf1.gamma_,
284  Y2*pf2.B_ - Y1*pf1.B_
285  );
286 }
287 
288 
289 // ************************************************************************* //
scalar Z(scalar p, scalar T) const
Return compression factor [].
dictionary dict
static autoPtr< adiabaticPerfectFluid > New(const dictionary &dict)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
adiabaticPerfectFluid(const Specie &sp, const scalar p0, const scalar rho0, const scalar gamma, const scalar B)
Construct from components.
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
scalar S(const scalar p, const scalar T) const
Return entropy [J/(kg K)].
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))
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
A class for handling words, derived from string.
Definition: word.H:59
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
const dimensionedScalar Pstd
Standard pressure.
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Adiabatic perfect fluid equation of state.
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52