incompressiblePerfectGasI.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-2019 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 
27 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Specie>
33 (
34  const Specie& sp, const scalar pRef
35 )
36 :
37  Specie(sp),
38  pRef_(pRef)
39 {}
40 
41 
42 template<class Specie>
44 (
45  const word& name,
47 )
48 :
49  Specie(name, ipg),
50  pRef_(ipg.pRef_)
51 {}
52 
53 
54 template<class Specie>
57 {
59  (
61  );
62 }
63 
64 
65 template<class Specie>
68 (
69  const dictionary& dict
70 )
71 {
73  (
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
81 template<class Specie>
83 (
84  scalar p,
85  scalar T
86 ) const
87 {
88  return pRef_/(this->R()*T);
89 }
90 
91 
92 template<class Specie>
94 (
95  scalar p,
96  scalar T
97 ) const
98 {
99  return p/this->rho(p, T);
100 }
101 
102 
103 template<class Specie>
105 (
106  scalar p,
107  scalar T
108 ) const
109 {
110  return 0;
111 }
112 
113 
114 template<class Specie>
116 (
117  scalar p,
118  scalar T
119 ) const
120 {
121  return 0;
122 }
123 
124 
125 template<class Specie>
127 (
128  scalar p,
129  scalar T
130 ) const
131 {
132  return 0;
133 }
134 
135 
136 template<class Specie>
138 (
139  scalar p,
140  scalar T
141 ) const
142 {
143  return 0;
144 }
145 
146 
147 template<class Specie>
149 (
150  scalar p,
151  scalar T
152 ) const
153 {
154  return 0;
155 }
156 
157 
158 template<class Specie>
160 (
161  scalar p,
162  scalar T
163 ) const
164 {
165  return 0;
166 }
167 
168 
169 template<class Specie>
171 (
172  scalar p,
173  scalar T
174 ) const
175 {
176  return 0;
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
181 
182 template<class Specie>
183 inline void Foam::incompressiblePerfectGas<Specie>::operator+=
184 (
186 )
187 {
188  scalar Y1 = this->Y();
189  Specie::operator+=(ipg);
190 
191  if (mag(this->Y()) > small)
192  {
193  Y1 /= this->Y();
194  const scalar Y2 = ipg.Y()/this->Y();
195 
196  pRef_ = Y1*pRef_ + Y2*ipg.pRef_;
197  }
198 }
199 
200 
201 template<class Specie>
203 {
204  Specie::operator*=(s);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
209 
210 template<class Specie>
211 inline Foam::incompressiblePerfectGas<Specie> Foam::operator+
212 (
215 )
216 {
217  Specie sp
218  (
219  static_cast<const Specie&>(ipg1)
220  + static_cast<const Specie&>(ipg2)
221  );
222 
223  if (mag(sp.Y()) < small)
224  {
226  (
227  sp,
228  ipg1.pRef_
229  );
230  }
231  else
232  {
233  const scalar Y1 = ipg1.Y()/sp.Y();
234  const scalar Y2 = ipg2.Y()/sp.Y();
235 
237  (
238  sp,
239  Y1*ipg1.pRef_ + Y2*ipg2.pRef_
240  );
241  }
242 }
243 
244 
245 template<class Specie>
246 inline Foam::incompressiblePerfectGas<Specie> Foam::operator*
247 (
248  const scalar s,
250 )
251 {
253  (
254  s*static_cast<const Specie&>(ipg),
255  ipg.pRef_
256  );
257 }
258 
259 
260 template<class Specie>
261 inline Foam::incompressiblePerfectGas<Specie> Foam::operator==
262 (
265 )
266 {
267  Specie sp
268  (
269  static_cast<const Specie&>(ipg1)
270  == static_cast<const Specie&>(ipg2)
271  );
272 
273  const scalar Y1 = ipg1.Y()/sp.Y();
274  const scalar Y2 = ipg2.Y()/sp.Y();
275 
277  (
278  sp,
279  Y2*ipg2.pRef_ - Y1*ipg1.pRef_
280  );
281 }
282 
283 
284 // ************************************************************************* //
dictionary dict
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Incompressible gas equation of state using a constant reference pressure in the perfect gas equation ...
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
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 Z(scalar p, scalar T) const
Return compression factor [].
A class for handling words, derived from string.
Definition: word.H:59
incompressiblePerfectGas(const Specie &sp, const scalar pRef)
Construct from components.
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].
const volScalarField & T
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
#define R(A, B, C, D, E, F, K, M)
PtrList< volScalarField > & Y
scalar psi(scalar p, scalar T) const
Return compressibility rho/p [s^2/m^2].
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
autoPtr< incompressiblePerfectGas > clone() const
Construct and return a clone.
static autoPtr< incompressiblePerfectGas > New(const dictionary &dict)