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-2021 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  // ***HGW This needs to be added, it is not 0
111  return 0;
112 }
113 
114 
115 template<class Specie>
116 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Cp
117 (
118  scalar p,
119  scalar T
120 ) const
121 {
122  return 0;
123 }
124 
125 
126 template<class Specie>
127 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::E
128 (
129  scalar p,
130  scalar T
131 ) const
132 {
133  // ***HGW This needs to be added, it is H - p/rho
134  return 0;
135 }
136 
137 
138 template<class Specie>
139 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Cv
140 (
141  scalar p,
142  scalar T
143 ) const
144 {
145  return 0;
146 }
147 
148 
149 template<class Specie>
150 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Sp
151 (
152  scalar p,
153  scalar T
154 ) const
155 {
156  scalar n = 1 - 1.0/gamma_;
157  return
158  -pow(p0_ + B_, 1.0/gamma_)*(pow((p + B_), n) - pow((Pstd + B_), n))
159  /(rho0_*T*n);
160 }
161 
162 
163 template<class Specie>
164 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Sv
165 (
166  scalar p,
167  scalar T
168 ) const
169 {
171  return 0;
172 }
173 
174 
175 template<class Specie>
177 (
178  scalar p,
179  scalar T
180 ) const
181 {
182  return
183  (rho0_/(gamma_*(p0_ + B_)))
184  *pow((p + B_)/(p0_ + B_), 1.0/gamma_ - 1.0);
185 }
186 
187 
188 template<class Specie>
189 inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Z
190 (
191  scalar p,
192  scalar T
193 ) const
194 {
195  return p/(rho(p, T)*this->R()*T);
196 }
197 
198 
199 template<class Specie>
201 (
202  scalar p,
203  scalar T
204 ) const
205 {
206  return 0;
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
211 
212 template<class Specie>
213 inline void Foam::adiabaticPerfectFluid<Specie>::operator+=
214 (
216 )
217 {
219 }
220 
221 
222 template<class Specie>
224 {
225  Specie::operator*=(s);
226 }
227 
228 
229 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
230 
231 template<class Specie>
232 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator+
233 (
236 )
237 {
239  return pf1;
240 }
241 
242 
243 template<class Specie>
244 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator*
245 (
246  const scalar s,
248 )
249 {
251  (
252  s*static_cast<const Specie&>(pf),
253  pf.p0_,
254  pf.rho0_,
255  pf.gamma_,
256  pf.B_
257  );
258 }
259 
260 
261 template<class Specie>
262 inline Foam::adiabaticPerfectFluid<Specie> Foam::operator==
263 (
266 )
267 {
269  return pf1;
270 }
271 
272 
273 // ************************************************************************* //
scalar Z(scalar p, scalar T) const
Return compression factor [].
dictionary dict
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
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:156
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 [s^2/m^2].
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 contribution [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].
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
const dimensionedScalar Pstd
Standard pressure.
#define noCoefficientMixing(Type)
Definition: specie.H:159
const volScalarField & T
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define R(A, B, C, D, E, F, K, M)
Adiabatic perfect fluid equation of state.
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
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
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370