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-2023 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  \verbatim
31  1/rho = C[0] + C[1]*T + C[2]*sqr(T) - C[3]*p - C[4]*p*T
32  \endverbatim
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 specification of the rPolynomial 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 (
79  const rPolynomial<Specie>&,
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 Classes
115 
116  //- Coefficient list class
117  class coeffList
118  :
119  public VectorSpace<coeffList, scalar, 5>
120  {
121  public:
122 
123  // Constructors
124 
125  //- Construct null
126  inline coeffList()
127  {}
128 
129  //- Construct from Istream
130  inline coeffList(Istream& is)
131  :
133  {}
134  };
135 
136 
137  // Private Data
138 
139  //- Density coefficients
140  coeffList C_;
141 
142 
143 public:
144 
145  // Constructors
146 
147  //- Construct from components
148  inline rPolynomial(const Specie& sp, const coeffList& coeffs);
149 
150  //- Construct from name and dictionary
151  rPolynomial(const word& name, const dictionary& dict);
152 
153  //- Construct as named copy
154  inline rPolynomial(const word& name, const rPolynomial&);
155 
156  //- Construct and return a clone
157  inline autoPtr<rPolynomial> clone() const;
158 
159 
160  // Member Functions
161 
162  //- Return the instantiated type name
163  static word typeName()
164  {
165  return "rPolynomial<" + word(Specie::typeName_()) + '>';
166  }
167 
168 
169  // Fundamental properties
170 
171  //- Is the equation of state is incompressible i.e. rho != f(p)
172  static const bool incompressible = false;
173 
174  //- Is the equation of state is isochoric i.e. rho = const
175  static const bool isochoric = false;
176 
177  //- Return density [kg/m^3]
178  inline scalar rho(scalar p, scalar T) const;
179 
180  //- Return enthalpy contribution [J/kg]
181  inline scalar H(const scalar p, const scalar T) const;
182 
183  //- Return Cp contribution [J/(kg K]
184  inline scalar Cp(scalar p, scalar T) const;
185 
186  //- Return internal energy contribution [J/kg]
187  inline scalar E(const scalar p, const scalar T) const;
188 
189  //- Return Cv contribution [J/(kg K]
190  inline scalar Cv(scalar p, scalar T) const;
191 
192  //- Return entropy contribution to the integral of Cp/T [J/kg/K]
193  inline scalar Sp(const scalar p, const scalar T) const;
194 
195  //- Return entropy contribution to the integral of Cv/T [J/kg/K]
196  inline scalar Sv(const scalar p, const scalar T) const;
197 
198  //- Return compressibility [s^2/m^2]
199  inline scalar psi(scalar p, scalar T) const;
200 
201  //- Return compression factor []
202  inline scalar Z(scalar p, scalar T) const;
203 
204  //- Return (Cp - Cv) [J/(kg K]
205  inline scalar CpMCv(scalar p, scalar T) const;
206 
207  //- Return volumetric coefficient of thermal expansion [1/T]
208  inline scalar alphav(const scalar p, const scalar T) const;
209 
210 
211  // IO
212 
213  //- Write to Ostream
214  void write(Ostream& os) const;
215 
216 
217  // Member Operators
218 
219  inline void operator+=(const rPolynomial&);
220  inline void operator*=(const scalar);
221 
222 
223  // Friend operators
224 
225  friend rPolynomial operator+ <Specie>
226  (
227  const rPolynomial&,
228  const rPolynomial&
229  );
230 
231  friend rPolynomial operator* <Specie>
232  (
233  const scalar s,
234  const rPolynomial&
235  );
236 
237  friend rPolynomial operator== <Specie>
238  (
239  const rPolynomial&,
240  const rPolynomial&
241  );
242 
243 
244  // Ostream Operator
245 
246  friend Ostream& operator<< <Specie>
247  (
248  Ostream&,
249  const rPolynomial&
250  );
251 };
252 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace Foam
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 #include "rPolynomialI.H"
261 
262 #ifdef NoRepository
263  #include "rPolynomial.C"
264 #endif
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #endif
269 
270 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Templated vector space.
Definition: VectorSpace.H:79
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Reciprocal polynomial equation of state for liquids and solids.
Definition: rPolynomial.H:118
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: rPolynomialI.H:93
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: rPolynomialI.H:107
void operator+=(const rPolynomial &)
Definition: rPolynomialI.H:146
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: rPolynomialI.H:86
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: rPolynomialI.H:115
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: rPolynomialI.H:72
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: rPolynomialI.H:136
static word typeName()
Return the instantiated type name.
Definition: rPolynomial.H:168
void write(Ostream &os) const
Write to Ostream.
Definition: rPolynomial.C:46
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: rPolynomialI.H:65
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: rPolynomialI.H:129
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: rPolynomialI.H:79
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: rPolynomialI.H:100
static const bool isochoric
Is the equation of state is isochoric i.e. rho = const.
Definition: rPolynomial.H:180
rPolynomial(const Specie &sp, const coeffList &coeffs)
Construct from components.
Definition: rPolynomialI.H:32
static const bool incompressible
Is the equation of state is incompressible i.e. rho != f(p)
Definition: rPolynomial.H:177
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: rPolynomialI.H:122
autoPtr< rPolynomial > clone() const
Construct and return a clone.
Definition: rPolynomialI.H:56
void operator*=(const scalar)
Definition: rPolynomialI.H:161
A class for handling words, derived from string.
Definition: word.H:62
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p