rPolynomial.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) 2019-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 Class
25  Foam::rPolynomial
26 
27 Description
28  Reciprocal polynomial equation of state for liquids and solids
29 
30  \f[
31  1/\rho = C_0 + C_1 T + C_2 T^2 - C_3 p - C_4 p T
32  \f]
33 
34  This polynomial for the reciprocal of the density provides a much better fit
35  than the equivalent polynomial for the density and has the advantage that it
36  support coefficient mixing to support liquid and solid mixtures in an
37  efficient manner.
38 
39 Usage
40  \table
41  Property | Description
42  C | Density polynomial coefficients
43  \endtable
44 
45  Example of the specification of the equation of state for pure water:
46  \verbatim
47  equationOfState
48  {
49  C (0.001278 -2.1055e-06 3.9689e-09 4.3772e-13 -2.0225e-16);
50  }
51  \endverbatim
52  Note: This fit is based on the small amount of data which is freely
53  available for the range 20-65degC and 1-100bar.
54 
55 SourceFiles
56  rPolynomialI.H
57  rPolynomial.C
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #ifndef rPolynomial_H
62 #define rPolynomial_H
63 
64 #include "autoPtr.H"
65 #include "VectorSpace.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of friend functions and operators
73 
74 template<class Specie> class rPolynomial;
75 
76 template<class Specie>
77 inline rPolynomial<Specie> operator+
78 (
80  const rPolynomial<Specie>&
81 );
82 
83 template<class Specie>
84 inline rPolynomial<Specie> operator*
85 (
86  const scalar,
87  const rPolynomial<Specie>&
88 );
89 
90 template<class Specie>
91 inline rPolynomial<Specie> operator==
92 (
93  const rPolynomial<Specie>&,
94  const rPolynomial<Specie>&
95 );
96 
97 template<class Specie>
98 Ostream& operator<<
99 (
100  Ostream&,
101  const rPolynomial<Specie>&
102 );
103 
104 
105 /*---------------------------------------------------------------------------*\
106  Class rPolynomial Declaration
107 \*---------------------------------------------------------------------------*/
108 
109 template<class Specie>
110 class rPolynomial
111 :
112  public Specie
113 {
114  // Private Data
115 
116  class coeffList
117  :
118  public VectorSpace<coeffList, scalar, 5>
119  {
120  public:
121 
122  // Constructors
123 
124  //- Construct null
125  inline coeffList()
126  {}
127 
128  //- Construct from Istream
129  inline coeffList(Istream& is)
130  :
132  {}
133  };
134 
135 
136  //- Density coefficients
137  coeffList C_;
138 
139 
140 public:
141 
142  // Constructors
143 
144  //- Construct from components
145  inline rPolynomial
146  (
147  const Specie& sp,
148  const coeffList& coeffs
149  );
150 
151  //- Construct from dictionary
152  rPolynomial(const dictionary& dict);
153 
154  //- Construct as named copy
155  inline rPolynomial(const word& name, const rPolynomial&);
156 
157  //- Construct and return a clone
158  inline autoPtr<rPolynomial> clone() const;
159 
160  // Selector from dictionary
161  inline static autoPtr<rPolynomial> New(const dictionary& dict);
162 
163 
164  // Member Functions
165 
166  //- Return the instantiated type name
167  static word typeName()
168  {
169  return "rPolynomial<" + word(Specie::typeName_()) + '>';
170  }
171 
173  // Fundamental properties
174 
175  //- Is the equation of state is incompressible i.e. rho != f(p)
176  static const bool incompressible = false;
177 
178  //- Is the equation of state is isochoric i.e. rho = const
179  static const bool isochoric = false;
180 
181  //- Return density [kg/m^3]
182  inline scalar rho(scalar p, scalar T) const;
183 
184  //- Return enthalpy contribution [J/kg]
185  inline scalar H(const scalar p, const scalar T) const;
186 
187  //- Return Cp contribution [J/(kg K]
188  inline scalar Cp(scalar p, scalar T) const;
189 
190  //- Return internal energy contribution [J/kg]
191  inline scalar E(const scalar p, const scalar T) const;
192 
193  //- Return Cv contribution [J/(kg K]
194  inline scalar Cv(scalar p, scalar T) const;
195 
196  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
197  inline scalar Sp(const scalar p, const scalar T) const;
198 
199  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
200  inline scalar Sv(const scalar p, const scalar T) const;
201 
202  //- Return compressibility [s^2/m^2]
203  inline scalar psi(scalar p, scalar T) const;
204 
205  //- Return compression factor []
206  inline scalar Z(scalar p, scalar T) const;
207 
208  //- Return (Cp - Cv) [J/(kg K]
209  inline scalar CpMCv(scalar p, scalar T) const;
210 
211 
212  // IO
213 
214  //- Write to Ostream
215  void write(Ostream& os) const;
216 
217 
218  // Member Operators
219 
220  inline void operator+=(const rPolynomial&);
221  inline void operator*=(const scalar);
222 
223 
224  // Friend operators
225 
226  friend rPolynomial operator+ <Specie>
227  (
228  const rPolynomial&,
229  const rPolynomial&
230  );
231 
232  friend rPolynomial operator* <Specie>
233  (
234  const scalar s,
235  const rPolynomial&
236  );
237 
238  friend rPolynomial operator== <Specie>
239  (
240  const rPolynomial&,
241  const rPolynomial&
242  );
243 
244 
245  // Ostream Operator
246 
247  friend Ostream& operator<< <Specie>
248  (
249  Ostream&,
250  const rPolynomial&
251  );
252 };
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace Foam
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #include "rPolynomialI.H"
262 
263 #ifdef NoRepository
264  #include "rPolynomial.C"
265 #endif
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #endif
270 
271 // ************************************************************************* //
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: rPolynomial.H:181
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: rPolynomialI.H:135
Reciprocal polynomial equation of state for liquids and solids.
Definition: rPolynomial.H:79
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Templated vector space.
Definition: VectorSpace.H:53
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: rPolynomialI.H:106
static word typeName()
Return the instantiated type name.
Definition: rPolynomial.H:172
autoPtr< rPolynomial > clone() const
Construct and return a clone.
Definition: rPolynomialI.H:58
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: rPolynomialI.H:85
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))
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: rPolynomial.H:184
A class for handling words, derived from string.
Definition: word.H:59
void operator+=(const rPolynomial &)
Definition: rPolynomialI.H:152
void operator*=(const scalar)
Definition: rPolynomialI.H:167
void write(Ostream &os) const
Write to Ostream.
Definition: rPolynomial.C:42
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: rPolynomialI.H:128
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: rPolynomialI.H:142
static autoPtr< rPolynomial > New(const dictionary &dict)
Definition: rPolynomialI.H:67
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: rPolynomialI.H:78
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 E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: rPolynomialI.H:99
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: rPolynomialI.H:120
volScalarField & p
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: rPolynomialI.H:92
rPolynomial(const Specie &sp, const coeffList &coeffs)
Construct from components.
Definition: rPolynomialI.H:32
Namespace for OpenFOAM.
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: rPolynomialI.H:113