perfectGasI.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) 2011-2020 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 "perfectGas.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
31 inline Foam::perfectGas<Specie>::perfectGas(const Specie& sp)
32 :
33  Specie(sp)
34 {}
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class Specie>
41 (
42  const word& name,
43  const perfectGas<Specie>& pg
44 )
45 :
46  Specie(name, pg)
47 {}
48 
49 
50 template<class Specie>
53 {
55 }
56 
57 
58 template<class Specie>
60 (
61  const dictionary& dict
62 )
63 {
65 }
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class Specie>
71 inline Foam::scalar Foam::perfectGas<Specie>::rho(scalar p, scalar T) const
72 {
73  return p/(this->R()*T);
74 }
75 
76 
77 template<class Specie>
78 inline Foam::scalar Foam::perfectGas<Specie>::H(scalar p, scalar T) const
79 {
80  return 0;
81 }
82 
83 
84 template<class Specie>
85 inline Foam::scalar Foam::perfectGas<Specie>::Cp(scalar p, scalar T) const
86 {
87  return 0;
88 }
89 
90 
91 template<class Specie>
92 inline Foam::scalar Foam::perfectGas<Specie>::E(scalar p, scalar T) const
93 {
94  return 0;
95 }
96 
97 
98 template<class Specie>
99 inline Foam::scalar Foam::perfectGas<Specie>::Cv(scalar p, scalar T) const
100 {
101  return 0;
102 }
103 
104 
105 template<class Specie>
106 inline Foam::scalar Foam::perfectGas<Specie>::Sp(scalar p, scalar T) const
107 {
108  return -this->R()*log(p/Pstd);
109 }
110 
111 
112 template<class Specie>
113 inline Foam::scalar Foam::perfectGas<Specie>::Sv(scalar p, scalar T) const
114 {
116  return 0;
117 }
118 
119 
120 template<class Specie>
121 inline Foam::scalar Foam::perfectGas<Specie>::psi(scalar p, scalar T) const
122 {
123  return 1.0/(this->R()*T);
124 }
125 
126 
127 template<class Specie>
128 inline Foam::scalar Foam::perfectGas<Specie>::Z(scalar p, scalar T) const
129 {
130  return 1;
131 }
132 
133 
134 template<class Specie>
135 inline Foam::scalar Foam::perfectGas<Specie>::CpMCv(scalar p, scalar T) const
136 {
137  return this->R();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
142 
143 template<class Specie>
145 {
146  Specie::operator+=(pg);
147 }
148 
149 
150 template<class Specie>
151 inline void Foam::perfectGas<Specie>::operator*=(const scalar s)
152 {
153  Specie::operator*=(s);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
158 
159 template<class Specie>
160 inline Foam::perfectGas<Specie> Foam::operator+
161 (
162  const perfectGas<Specie>& pg1,
163  const perfectGas<Specie>& pg2
164 )
165 {
166  return perfectGas<Specie>
167  (
168  static_cast<const Specie&>(pg1) + static_cast<const Specie&>(pg2)
169  );
170 }
171 
172 
173 template<class Specie>
174 inline Foam::perfectGas<Specie> Foam::operator*
175 (
176  const scalar s,
177  const perfectGas<Specie>& pg
178 )
179 {
180  return perfectGas<Specie>(s*static_cast<const Specie&>(pg));
181 }
182 
183 
184 template<class Specie>
185 inline Foam::perfectGas<Specie> Foam::operator==
186 (
187  const perfectGas<Specie>& pg1,
188  const perfectGas<Specie>& pg2
189 )
190 {
191  return perfectGas<Specie>
192  (
193  static_cast<const Specie&>(pg1) == static_cast<const Specie&>(pg2)
194  );
195 }
196 
197 
198 // ************************************************************************* //
void operator*=(const scalar)
Definition: perfectGasI.H:151
dictionary dict
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: perfectGasI.H:135
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
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: perfectGasI.H:113
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: perfectGasI.H:106
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: perfectGasI.H:99
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: perfectGasI.H:78
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: perfectGasI.H:85
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: perfectGasI.H:121
void operator+=(const perfectGas &)
Definition: perfectGasI.H:144
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))
autoPtr< perfectGas > clone() const
Construct and return a clone.
Definition: perfectGasI.H:52
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].
Definition: perfectGasI.H:92
Perfect gas equation of state.
Definition: perfectGas.H:47
const dimensionedScalar Pstd
Standard pressure.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static autoPtr< perfectGas > New(const dictionary &dict)
Definition: perfectGasI.H:60
#define R(A, B, C, D, E, F, K, M)
perfectGas(const Specie &sp)
Construct from components.
Definition: perfectGasI.H:31
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: perfectGasI.H:128
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: perfectGasI.H:71