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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31
32template<class Specie>
34(
35 const Specie& sp,
36 const scalar p0,
37 const scalar rho0,
38 const scalar gamma,
39 const scalar B
40)
41:
42 Specie(sp),
43 p0_(p0),
44 rho0_(rho0),
45 gamma_(gamma),
46 B_(B)
47{}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52template<class Specie>
54(
55 const word& name,
57)
58:
59 Specie(name, pf),
60 p0_(pf.p0_),
61 rho0_(pf.rho0_),
62 gamma_(pf.gamma_),
63 B_(pf.B_)
64{}
65
66
67template<class Specie>
70{
72}
73
74
75template<class Specie>
78(
79 const dictionary& dict
80)
81{
83}
84
85
86// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87
88template<class Specie>
90(
91 scalar p,
92 scalar T
93) const
94{
95 return rho0_*pow((p + B_)/(p0_ + B_), 1.0/gamma_);
96}
97
98
99template<class Specie>
101(
102 scalar p,
103 scalar T
104) const
105{
106 return 0;
107}
108
109
110template<class Specie>
112(
113 scalar p,
114 scalar T
115) const
116{
117 return 0;
118}
119
120
121template<class Specie>
123(
124 scalar p,
125 scalar T
126) const
127{
128 return 0;
129}
130
131
132template<class Specie>
134(
135 scalar p,
136 scalar T
137) const
138{
139 return 0;
140}
141
142
143template<class Specie>
145(
146 scalar p,
147 scalar T
148) const
149{
150 scalar n = 1 - 1.0/gamma_;
151 return
152 -pow(p0_ + B_, 1.0/gamma_)*(pow((p + B_), n) - pow((Pstd + B_), n))
153 /(rho0_*T*n);
154}
155
156
157template<class Specie>
159(
160 scalar p,
161 scalar T
162) const
163{
164 return
165 (rho0_/(gamma_*(p0_ + B_)))
166 *pow((p + B_)/(p0_ + B_), 1.0/gamma_ - 1.0);
167}
168
169
170template<class Specie>
171inline Foam::scalar Foam::adiabaticPerfectFluid<Specie>::Z(scalar, scalar) const
172{
173 return 1;
174}
175
176
177template<class Specie>
179(
180 scalar p,
181 scalar T
182) const
183{
184 return 0;
185}
186
187
188// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
189
190template<class Specie>
192(
194)
195{
196 scalar Y1 = this->Y();
197 Specie::operator+=(pf);
198
199 if (mag(this->Y()) > SMALL)
200 {
201 Y1 /= this->Y();
202 const scalar Y2 = pf.Y()/this->Y();
203
204 p0_ = Y1*p0_ + Y2*pf.p0_;
205 rho0_ = Y1*rho0_ + Y2*pf.rho0_;
206 gamma_ = Y1*gamma_ + Y2*pf.gamma_;
207 B_ = Y1*B_ + Y2*pf.B_;
208 }
209}
210
211
212template<class Specie>
214{
215 Specie::operator*=(s);
216}
217
218
219// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
220
221template<class Specie>
222inline Foam::adiabaticPerfectFluid<Specie> Foam::operator+
223(
226)
227{
228 Specie sp
229 (
230 static_cast<const Specie&>(pf1)
231 + static_cast<const Specie&>(pf2)
232 );
233
234 if (mag(sp.Y()) < SMALL)
235 {
237 (
238 sp,
239 pf1.p0_,
240 pf1.rho0_,
241 pf1.gamma_,
242 pf1.B_
243 );
244 }
245 else
246 {
247 const scalar Y1 = pf1.Y()/sp.Y();
248 const scalar Y2 = pf2.Y()/sp.Y();
249
250 return adiabaticPerfectFluid<Specie>
251 (
252 sp,
253 Y1*pf1.p0_ + Y2*pf2.p0_,
254 Y1*pf1.rho0_ + Y2*pf2.rho0_,
255 Y1*pf1.gamma_ + Y2*pf2.gamma_,
256 Y1*pf1.B_ + Y2*pf2.B_
257 );
258 }
259}
260
261
262template<class Specie>
263inline Foam::adiabaticPerfectFluid<Specie> Foam::operator*
264(
265 const scalar s,
266 const adiabaticPerfectFluid<Specie>& pf
267)
268{
269 return adiabaticPerfectFluid<Specie>
270 (
271 s*static_cast<const Specie&>(pf),
272 pf.p0_,
273 pf.rho0_,
274 pf.gamma_,
275 pf.B_
276 );
277}
278
279
280template<class Specie>
281inline Foam::adiabaticPerfectFluid<Specie> Foam::operator==
282(
283 const adiabaticPerfectFluid<Specie>& pf1,
284 const adiabaticPerfectFluid<Specie>& pf2
285)
286{
287 Specie sp
288 (
289 static_cast<const Specie&>(pf1)
290 == static_cast<const Specie&>(pf2)
291 );
292
293 const scalar Y1 = pf1.Y()/sp.Y();
294 const scalar Y2 = pf2.Y()/sp.Y();
295
296 return adiabaticPerfectFluid<Specie>
297 (
298 sp,
299 Y2*pf2.p0_ - Y1*pf1.p0_,
300 Y2*pf2.rho0_ - Y1*pf1.rho0_,
301 Y2*pf2.gamma_ - Y1*pf1.gamma_,
302 Y2*pf2.B_ - Y1*pf1.B_
303 );
304}
305
306
307// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
label n
compactSpatialTensor S
The joint motion sub-space (3-DoF)
Definition: joint.H:132
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Adiabatic perfect fluid equation of state.
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
autoPtr< adiabaticPerfectFluid > clone() const
Construct and return a clone.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & psi
const volScalarField & T
const volScalarField & Cv
Definition: EEqn.H:8
const scalar gamma
Definition: EEqn.H:9
const volScalarField & Cp
Definition: EEqn.H:7
const volScalarField & p0
Definition: EEqn.H:36
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
scalar rho0
dictionary dict