Loading...
Searching...
No Matches
rPolynomialI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "rPolynomial.H"
29
30// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31
32template<class Specie>
34(
35 const Specie& sp,
36 const coeffList& coeffs
37)
38:
39 Specie(sp),
40 C_(coeffs)
41{}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46template<class Specie>
48(
49 const word& name,
51)
52:
53 Specie(name, rp),
54 C_(rp.C_)
55{}
56
57
58template<class Specie>
63}
64
65
66template<class Specie>
69(
70 const dictionary& dict
71)
72{
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
79template<class Specie>
80inline Foam::scalar Foam::rPolynomial<Specie>::rho(scalar p, scalar T) const
81{
82 return 1/(C_[0] + (C_[1] + C_[2]*T - C_[4]*p)*T - C_[3]*p);
83}
84
85
86template<class Specie>
87inline Foam::scalar Foam::rPolynomial<Specie>::H(scalar p, scalar T) const
88{
89 return 0;
90}
91
92
93template<class Specie>
94inline Foam::scalar Foam::rPolynomial<Specie>::Cp(scalar p, scalar T) const
95{
96 return 0;
97}
98
99
100template<class Specie>
101inline Foam::scalar Foam::rPolynomial<Specie>::E(scalar p, scalar T) const
102{
103 return 0;
104}
105
106
107template<class Specie>
108inline Foam::scalar Foam::rPolynomial<Specie>::Cv(scalar p, scalar T) const
109{
110 return 0;
111}
112
113
114template<class Specie>
115inline Foam::scalar Foam::rPolynomial<Specie>::S(scalar p, scalar T) const
116{
117 return 0;
118}
119
120
121template<class Specie>
122inline Foam::scalar Foam::rPolynomial<Specie>::psi(scalar p, scalar T) const
123{
124 return sqr(rho(p, T))*(C_[3] + C_[4]*T);
125}
126
127
128template<class Specie>
129inline Foam::scalar Foam::rPolynomial<Specie>::Z(scalar p, scalar T) const
130{
131 return 1;
132}
133
134
135template<class Specie>
136inline Foam::scalar Foam::rPolynomial<Specie>::CpMCv(scalar p, scalar T) const
137{
138 return 0;
139}
140
141
142// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
143
144template<class Specie>
146(
148)
149{
150 const scalar Y1 = this->Y();
151 Specie::operator+=(rp);
152
153 if (mag(this->Y()) > SMALL)
155 C_ = (Y1*C_ + rp.Y()*rp.C_)/this->Y();
156 }
157}
158
159
160template<class Specie>
161inline void Foam::rPolynomial<Specie>::operator*=(const scalar s)
162{
163 Specie::operator*=(s);
164}
165
166
167// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
168
169template<class Specie>
170inline Foam::rPolynomial<Specie> Foam::operator+
171(
172 const rPolynomial<Specie>& rp1,
173 const rPolynomial<Specie>& rp2
174)
175{
176 Specie sp
177 (
178 static_cast<const Specie&>(rp1)
179 + static_cast<const Specie&>(rp2)
180 );
181
182 if (mag(sp.Y()) < SMALL)
183 {
185 (
186 sp,
187 rp1.C_
188 );
189 }
190
192 (
193 sp,
194 (rp1.Y()*rp1.C_ + rp2.Y()*rp2.C_)/sp.Y()
195 );
196}
197
198
199template<class Specie>
200inline Foam::rPolynomial<Specie> Foam::operator*
201(
202 const scalar s,
204)
205{
207 (
208 s*static_cast<const Specie&>(rp),
209 rp.C_
210 );
211}
212
213
214template<class Specie>
215inline Foam::rPolynomial<Specie> Foam::operator==
216(
217 const rPolynomial<Specie>& rp1,
218 const rPolynomial<Specie>& rp2
219)
220{
221 Specie sp
222 (
223 static_cast<const Specie&>(rp1)
224 == static_cast<const Specie&>(rp2)
225 );
226
228 (
229 sp,
230 (rp1.Y()*rp1.C_ - rp2.Y()*rp2.C_)/sp.Y()
231 );
232}
233
234
235// ************************************************************************* //
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Reciprocal polynomial equation of state for liquids and solids.
scalar Cv(scalar p, scalar T) const
Return Cv departure [J/(kg K].
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
static autoPtr< rPolynomial > New(const dictionary &dict)
rPolynomial(const Specie &sp, const coeffList &coeffs)
Construct from components.
scalar Z(scalar p, scalar T) const
Return compression factor [].
autoPtr< rPolynomial > clone() const
Construct and return a clone.
void operator*=(const scalar)
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
PtrList< volScalarField > & Y
regionProperties rp(runTime)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict