perfectFluidI.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) 2012-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 "perfectFluid.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const scalar R,
35  const scalar rho0
36 )
37 :
38  Specie(sp),
39  R_(R),
40  rho0_(rho0)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Specie>
48 (
49  const word& name,
50  const perfectFluid<Specie>& pf
51 )
52 :
53  Specie(name, pf),
54  R_(pf.R_),
55  rho0_(pf.rho0_)
56 {}
57 
58 
59 template<class Specie>
62 {
64 }
65 
66 
67 template<class Specie>
70 (
71  const dictionary& dict
72 )
73 {
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class Specie>
81 inline Foam::scalar Foam::perfectFluid<Specie>::R() const
82 {
83  return R_;
84 }
85 
86 
87 template<class Specie>
88 inline Foam::scalar Foam::perfectFluid<Specie>::rho(scalar p, scalar T) const
89 {
90  return rho0_ + p/(this->R()*T);
91 }
92 
93 
94 template<class Specie>
95 inline Foam::scalar Foam::perfectFluid<Specie>::H(scalar p, scalar T) const
96 {
97  return p/rho(p, T) - Pstd/rho(Pstd, T);
98 }
99 
100 
101 template<class Specie>
102 inline Foam::scalar Foam::perfectFluid<Specie>::Cp(scalar p, scalar T) const
103 {
104  return
105  sqr(p/(rho(p, T)*T))/this->R()
106  - sqr(Pstd/(rho(Pstd, T)*T))/this->R();
107 }
108 
109 
110 template<class Specie>
111 inline Foam::scalar Foam::perfectFluid<Specie>::E(scalar p, scalar T) const
112 {
113  return 0;
114 }
115 
116 
117 template<class Specie>
118 inline Foam::scalar Foam::perfectFluid<Specie>::Cv(scalar p, scalar T) const
119 {
120  return 0;
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::perfectFluid<Specie>::Sp(scalar p, scalar T) const
126 {
127  return -this->R()*log(p/Pstd);
128 }
129 
130 
131 template<class Specie>
132 inline Foam::scalar Foam::perfectFluid<Specie>::Sv(scalar p, scalar T) const
133 {
135  return 0;
136 }
137 
138 
139 template<class Specie>
140 inline Foam::scalar Foam::perfectFluid<Specie>::psi(scalar p, scalar T) const
141 {
142  return 1.0/(this->R()*T);
143 }
144 
145 
146 template<class Specie>
147 inline Foam::scalar Foam::perfectFluid<Specie>::Z(scalar p, scalar T) const
148 {
149  return p/(rho(p, T)*this->R()*T);
150 }
151 
152 
153 template<class Specie>
154 inline Foam::scalar Foam::perfectFluid<Specie>::CpMCv(scalar p, scalar T) const
155 {
156  const scalar R = this->R();
157  const scalar rho = this->rho(p, T);
158 
159  return R*sqr(p/(rho*R*T));
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
164 
165 template<class Specie>
166 inline void Foam::perfectFluid<Specie>::operator+=
167 (
168  const perfectFluid<Specie>& pf
169 )
170 {
172 }
173 
174 
175 template<class Specie>
176 inline void Foam::perfectFluid<Specie>::operator*=(const scalar s)
177 {
178  Specie::operator*=(s);
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
183 
184 template<class Specie>
185 inline Foam::perfectFluid<Specie> Foam::operator+
186 (
187  const perfectFluid<Specie>& pf1,
188  const perfectFluid<Specie>& pf2
189 )
190 {
192  return pf1;
193 }
194 
195 
196 template<class Specie>
197 inline Foam::perfectFluid<Specie> Foam::operator*
198 (
199  const scalar s,
200  const perfectFluid<Specie>& pf
201 )
202 {
203  return perfectFluid<Specie>
204  (
205  s*static_cast<const Specie&>(pf),
206  pf.R_,
207  pf.rho0_
208  );
209 }
210 
211 
212 template<class Specie>
213 inline Foam::perfectFluid<Specie> Foam::operator==
214 (
215  const perfectFluid<Specie>& pf1,
216  const perfectFluid<Specie>& pf2
217 )
218 {
220  return pf1;
221 }
222 
223 
224 // ************************************************************************* //
autoPtr< perfectFluid > clone() const
Construct and return a clone.
Definition: perfectFluidI.H:61
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedSymmTensor sqr(const dimensionedVector &dv)
perfectFluid(const Specie &sp, const scalar R, const scalar rho0)
Construct from components.
Definition: perfectFluidI.H:32
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectFluidI.H:88
static autoPtr< perfectFluid > New(const dictionary &dict)
Definition: perfectFluidI.H:70
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:81
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))
A class for handling words, derived from string.
Definition: word.H:59
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
void operator*=(const scalar)
scalar Z(scalar p, scalar T) const
Return compression factor [].
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: perfectFluidI.H:95
const dimensionedScalar Pstd
Standard pressure.
#define noCoefficientMixing(Type)
Definition: specie.H:159
const volScalarField & T
Simple extension of the perfect gas equation of state to liquids by the addition of a constant densit...
Definition: perfectFluid.H:54
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
#define R(A, B, C, D, E, F, K, M)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].