rPolynomialI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "rPolynomial.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const coeffList& coeffs
35 )
36 :
37  Specie(sp),
38  C_(coeffs)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class Specie>
46 (
47  const word& name,
48  const rPolynomial<Specie>& rp
49 )
50 :
51  Specie(name, rp),
52  C_(rp.C_)
53 {}
54 
55 
56 template<class Specie>
59 {
61 }
62 
63 
64 template<class Specie>
67 (
68  const dictionary& dict
69 )
70 {
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
77 template<class Specie>
78 inline Foam::scalar Foam::rPolynomial<Specie>::rho(scalar p, scalar T) const
79 {
80  return 1/(C_[0] + (C_[1] + C_[2]*T - C_[4]*p)*T - C_[3]*p);
81 }
82 
83 
84 template<class Specie>
85 inline Foam::scalar Foam::rPolynomial<Specie>::H(scalar p, scalar T) const
86 {
87  return 0;
88 }
89 
90 
91 template<class Specie>
92 inline Foam::scalar Foam::rPolynomial<Specie>::Cp(scalar p, scalar T) const
93 {
94  return 0;
95 }
96 
97 
98 template<class Specie>
99 inline Foam::scalar Foam::rPolynomial<Specie>::E(scalar p, scalar T) const
100 {
101  return 0;
102 }
103 
104 
105 template<class Specie>
106 inline Foam::scalar Foam::rPolynomial<Specie>::Cv(scalar p, scalar T) const
107 {
108  return 0;
109 }
110 
111 
112 template<class Specie>
113 inline Foam::scalar Foam::rPolynomial<Specie>::Sp(scalar p, scalar T) const
114 {
115  return 0;
116 }
117 
118 
119 template<class Specie>
120 inline Foam::scalar Foam::rPolynomial<Specie>::Sv(scalar p, scalar T) const
121 {
123  return 0;
124 }
125 
126 
127 template<class Specie>
128 inline Foam::scalar Foam::rPolynomial<Specie>::psi(scalar p, scalar T) const
129 {
130  return sqr(rho(p, T))*(C_[3] + C_[4]*T);
131 }
132 
133 
134 template<class Specie>
135 inline Foam::scalar Foam::rPolynomial<Specie>::Z(scalar p, scalar T) const
136 {
137  return 1;
138 }
139 
140 
141 template<class Specie>
142 inline Foam::scalar Foam::rPolynomial<Specie>::CpMCv(scalar p, scalar T) const
143 {
144  return 0;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
149 
150 template<class Specie>
151 inline void Foam::rPolynomial<Specie>::operator+=
152 (
153  const rPolynomial<Specie>& rp
154 )
155 {
156  const scalar Y1 = this->Y();
157  Specie::operator+=(rp);
158 
159  if (mag(this->Y()) > small)
160  {
161  C_ = (Y1*C_ + rp.Y()*rp.C_)/this->Y();
162  }
163 }
164 
165 
166 template<class Specie>
167 inline void Foam::rPolynomial<Specie>::operator*=(const scalar s)
168 {
169  Specie::operator*=(s);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
174 
175 template<class Specie>
176 inline Foam::rPolynomial<Specie> Foam::operator+
177 (
178  const rPolynomial<Specie>& rp1,
179  const rPolynomial<Specie>& rp2
180 )
181 {
182  Specie sp
183  (
184  static_cast<const Specie&>(rp1)
185  + static_cast<const Specie&>(rp2)
186  );
187 
188  if (mag(sp.Y()) < small)
189  {
190  return rPolynomial<Specie>
191  (
192  sp,
193  rp1.C_
194  );
195  }
196  else
197  {
198  return rPolynomial<Specie>
199  (
200  sp,
201  (rp1.Y()*rp1.C_ + rp2.Y()*rp2.C_)/sp.Y()
202  );
203  }
204 
205  return rp1;
206 }
207 
208 
209 template<class Specie>
210 inline Foam::rPolynomial<Specie> Foam::operator*
211 (
212  const scalar s,
213  const rPolynomial<Specie>& rp
214 )
215 {
216  return rPolynomial<Specie>
217  (
218  s*static_cast<const Specie&>(rp),
219  rp.C_
220  );
221 }
222 
223 
224 template<class Specie>
225 inline Foam::rPolynomial<Specie> Foam::operator==
226 (
227  const rPolynomial<Specie>& rp1,
228  const rPolynomial<Specie>& rp2
229 )
230 {
231  return rPolynomial<Specie>
232  (
233  static_cast<const Specie&>(rp1) == static_cast<const Specie&>(rp2),
235  );
236 }
237 
238 
239 // ************************************************************************* //
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Reciprocal polynomial equation of state for liquids and solids.
Definition: rPolynomial.H:79
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: rPolynomialI.H:106
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))
A class for handling words, derived from string.
Definition: word.H:59
regionProperties rp(runTime)
void operator*=(const scalar)
Definition: rPolynomialI.H:167
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: rPolynomialI.H:128
const volScalarField & T
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: rPolynomialI.H:142
PtrList< volScalarField > & Y
static autoPtr< rPolynomial > New(const dictionary &dict)
Definition: rPolynomialI.H:67
dimensioned< scalar > mag(const dimensioned< Type > &)
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
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
rPolynomial(const Specie &sp, const coeffList &coeffs)
Construct from components.
Definition: rPolynomialI.H:32
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