perfectFluidI.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) 2012-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
28#include "perfectFluid.H"
29#include "specie.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class Specie>
35(
36 const Specie& sp,
37 const scalar R,
38 const scalar rho0
39)
40:
41 Specie(sp),
42 R_(R),
43 rho0_(rho0)
44{}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49template<class Specie>
51(
52 const word& name,
53 const perfectFluid<Specie>& pf
54)
55:
56 Specie(name, pf),
57 R_(pf.R_),
58 rho0_(pf.rho0_)
59{}
60
61
62template<class Specie>
65{
67}
68
69
70template<class Specie>
73(
74 const dictionary& dict
75)
76{
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
83template<class Specie>
84inline Foam::scalar Foam::perfectFluid<Specie>::R() const
85{
86 return R_;
87}
88
89
90template<class Specie>
91inline Foam::scalar Foam::perfectFluid<Specie>::rho(scalar p, scalar T) const
92{
93 return rho0_ + p/(this->R()*T);
94}
95
96
97template<class Specie>
98inline Foam::scalar Foam::perfectFluid<Specie>::H(scalar p, scalar T) const
99{
100 return 0;
101}
102
103
104template<class Specie>
105inline Foam::scalar Foam::perfectFluid<Specie>::Cp(scalar p, scalar T) const
106{
107 return 0;
108}
109
110
111template<class Specie>
112inline Foam::scalar Foam::perfectFluid<Specie>::E(scalar p, scalar T) const
113{
114 return 0;
115}
116
117
118template<class Specie>
119inline Foam::scalar Foam::perfectFluid<Specie>::Cv(scalar p, scalar T) const
120{
121 return 0;
122}
123
124
125template<class Specie>
126inline Foam::scalar Foam::perfectFluid<Specie>::S(scalar p, scalar T) const
127{
128 return -this->R()*log(p/Pstd);
129}
130
131
132template<class Specie>
133inline Foam::scalar Foam::perfectFluid<Specie>::psi(scalar p, scalar T) const
134{
135 return 1.0/(this->R()*T);
136}
137
138
139template<class Specie>
140inline Foam::scalar Foam::perfectFluid<Specie>::Z(scalar p, scalar T) const
141{
142 return 1;
143}
144
145
146template<class Specie>
147inline Foam::scalar Foam::perfectFluid<Specie>::CpMCv(scalar p, scalar T) const
148{
149 const scalar R = this->R();
150 const scalar rho = this->rho(p, T);
151
152 return R*sqr(p/(rho*R*T));
153}
154
155
156// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157
158template<class Specie>
160(
161 const perfectFluid<Specie>& pf
162)
163{
164 scalar Y1 = this->Y();
165 Specie::operator+=(pf);
166
167 if (mag(this->Y()) > SMALL)
168 {
169 Y1 /= this->Y();
170 const scalar Y2 = pf.Y()/this->Y();
171
172 R_ = 1.0/(Y1/R_ + Y2/pf.R_);
173 rho0_ = Y1*rho0_ + Y2*pf.rho0_;
174 }
175}
176
177
178template<class Specie>
180{
181 Specie::operator*=(s);
182}
183
184
185// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
186
187template<class Specie>
188inline Foam::perfectFluid<Specie> Foam::operator+
189(
190 const perfectFluid<Specie>& pf1,
191 const perfectFluid<Specie>& pf2
192)
193{
194 Specie sp
195 (
196 static_cast<const Specie&>(pf1)
197 + static_cast<const Specie&>(pf2)
198 );
199
200 if (mag(sp.Y()) < SMALL)
201 {
203 (
204 sp,
205 pf1.R_,
206 pf1.rho0_
207 );
208 }
209 else
210 {
211 const scalar Y1 = pf1.Y()/sp.Y();
212 const scalar Y2 = pf2.Y()/sp.Y();
213
214 return perfectFluid<Specie>
215 (
216 sp,
217 1.0/(Y1/pf1.R_ + Y2/pf2.R_),
218 Y1*pf1.rho0_ + Y2*pf2.rho0_
219 );
220 }
221}
222
223
224template<class Specie>
225inline Foam::perfectFluid<Specie> Foam::operator*
226(
227 const scalar s,
228 const perfectFluid<Specie>& pf
229)
230{
231 return perfectFluid<Specie>
232 (
233 s*static_cast<const Specie&>(pf),
234 pf.R_,
235 pf.rho0_
236 );
237}
238
239
240template<class Specie>
241inline Foam::perfectFluid<Specie> Foam::operator==
242(
243 const perfectFluid<Specie>& pf1,
244 const perfectFluid<Specie>& pf2
245)
246{
247 Specie sp
248 (
249 static_cast<const Specie&>(pf1)
250 == static_cast<const Specie&>(pf2)
251 );
252
253 const scalar Y1 = pf1.Y()/sp.Y();
254 const scalar Y2 = pf2.Y()/sp.Y();
255
256 return perfectFluid<Specie>
257 (
258 sp,
259 1.0/(Y2/pf2.R_ - Y1/pf1.R_),
260 Y2*pf2.rho0_ - Y1*pf1.rho0_
261 );
262}
263
264
265// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
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
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
Perfect gas equation of state.
Definition: perfectFluid.H:91
autoPtr< perfectFluid > clone() const
Construct and return a clone.
Definition: perfectFluidI.H:64
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar R() const
Return fluid constant [J/(kg K)].
Definition: perfectFluidI.H:84
void operator*=(const scalar)
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 volScalarField & Cp
Definition: EEqn.H:7
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))
const scalar Pstd
Standard pressure: default in [Pa].
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log(const dimensionedScalar &ds)
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