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-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
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 template<class Specie>
50 (
51  const word& name,
53 )
54 :
55  Specie(name, pf),
56  p0_(pf.p0_),
57  rho0_(pf.rho0_),
58  gamma_(pf.gamma_),
59  B_(pf.B_)
60 {}
61 
62 
63 template<class Specie>
66 {
68  (
70  );
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
76 template<class Specie>
78 (
79  scalar p,
80  scalar T
81 ) const
82 {
83  return rho0_*pow((p + B_)/(p0_ + B_), 1.0/gamma_);
84 }
85 
86 
87 template<class Specie>
89 (
90  scalar p,
91  scalar T
92 ) const
93 {
94  // ***HGW This needs to be added, it is not 0
95  return 0;
96 }
97 
98 
99 template<class Specie>
101 (
102  scalar p,
103  scalar T
104 ) const
105 {
106  return 0;
107 }
108 
109 
110 template<class Specie>
112 (
113  scalar p,
114  scalar T
115 ) const
116 {
117  // ***HGW This needs to be added, it is H - p/rho
118  return 0;
119 }
120 
121 
122 template<class Specie>
124 (
125  scalar p,
126  scalar T
127 ) const
128 {
129  return 0;
130 }
131 
132 
133 template<class Specie>
135 (
136  scalar p,
137  scalar T
138 ) const
139 {
140  scalar n = 1 - 1.0/gamma_;
141  return
142  -pow(p0_ + B_, 1.0/gamma_)*(pow((p + B_), n) - pow((Pstd + B_), n))
143  /(rho0_*T*n);
144 }
145 
146 
147 template<class Specie>
149 (
150  scalar p,
151  scalar T
152 ) const
153 {
155  return 0;
156 }
157 
158 
159 template<class Specie>
161 (
162  scalar p,
163  scalar T
164 ) const
165 {
166  return
167  (rho0_/(gamma_*(p0_ + B_)))
168  *pow((p + B_)/(p0_ + B_), 1.0/gamma_ - 1.0);
169 }
170 
171 
172 template<class Specie>
174 (
175  scalar p,
176  scalar T
177 ) const
178 {
179  return p/(rho(p, T)*this->R()*T);
180 }
181 
182 
183 template<class Specie>
185 (
186  scalar p,
187  scalar T
188 ) const
189 {
190  return 0;
191 }
192 
193 
194 template<class Specie>
196 (
197  scalar p,
198  scalar T
199 ) const
200 {
201  return 0;
202 }
203 
204 
205 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
206 
207 template<class Specie>
209 (
211 )
212 {
214 }
215 
216 
217 template<class Specie>
219 {
220  Specie::operator*=(s);
221 }
222 
223 
224 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
225 
226 template<class Specie>
227 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator+
228 (
231 )
232 {
234  return pf1;
235 }
236 
237 
238 template<class Specie>
239 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator*
240 (
241  const scalar s,
242  const adiabaticPerfectFluid<Specie>& pf
243 )
244 {
245  return adiabaticPerfectFluid<Specie>
246  (
247  s*static_cast<const Specie&>(pf),
248  pf.p0_,
249  pf.rho0_,
250  pf.gamma_,
251  pf.B_
252  );
253 }
254 
255 
256 template<class Specie>
257 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator==
258 (
259  const adiabaticPerfectFluid<Specie>& pf1,
260  const adiabaticPerfectFluid<Specie>& pf2
261 )
262 {
263  noCoefficientMixing(adiabaticPerfectFluid);
264  return pf1;
265 }
266 
267 
268 // ************************************************************************* //
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
label n
Adiabatic perfect fluid equation of state for liquids:
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
adiabaticPerfectFluid(const Specie &sp, const scalar p0, const scalar rho0, const scalar gamma, const scalar B)
Construct from components.
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
scalar Z(scalar p, scalar T) const
Return compression factor [].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar Pstd
Standard pressure.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho0
volScalarField & p
#define noCoefficientMixing(Type)
Definition: specie.H:159