icoPolynomialI.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) 2011-2017 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 "icoPolynomial.H"
29
30// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31
32template<class Specie, int PolySize>
34(
35 const Specie& sp,
36 const Polynomial<PolySize>& rhoCoeffs
37)
38:
39 Specie(sp),
40 rhoCoeffs_(rhoCoeffs)
41{}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46template<class Specie, int PolySize>
48(
49 const word& name,
51)
52:
53 Specie(name, ip),
54 rhoCoeffs_(ip.rhoCoeffs_)
55{}
56
57
58template<class Specie, int PolySize>
61{
63}
64
65
66template<class Specie, int PolySize>
69{
71}
72
73
74// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75
76template<class Specie, int PolySize>
78(
79 scalar p,
80 scalar T
81) const
82{
83 return rhoCoeffs_.value(T);
84}
85
86
87template<class Specie, int PolySize>
89(
90 scalar p,
91 scalar T
92) const
93{
94 return 0;
95}
96
97
98template<class Specie, int PolySize>
100(
101 scalar p,
102 scalar T
103) const
104{
105 return 0;
106}
107
108
109template<class Specie, int PolySize>
111(
112 scalar p,
113 scalar T
114) const
115{
116 return 0;
117}
118
119
120template<class Specie, int PolySize>
122(
123 scalar p,
124 scalar T
125) const
126{
127 return 0;
128}
129
130
131template<class Specie, int PolySize>
133(
134 scalar p,
135 scalar T
136) const
137{
138 return 0;
139}
140
141
142template<class Specie, int PolySize>
144(
145 scalar p,
146 scalar T
147) const
148{
149 return 0;
150}
151
152
153template<class Specie, int PolySize>
155(
156 scalar p,
157 scalar T
158) const
159{
160 return 0;
161}
162
163
164template<class Specie, int PolySize>
166(
167 scalar p,
168 scalar T
169) const
170{
171 return 0;
172}
173
174
175// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
176
177template<class Specie, int PolySize>
179(
181)
182{
183 scalar Y1 = this->Y();
184 Specie::operator+=(ip);
185
186 if (mag(this->Y()) > SMALL)
187 {
188 Y1 /= this->Y();
189 const scalar Y2 = ip.Y()/this->Y();
190
191 rhoCoeffs_ = Y1*rhoCoeffs_ + Y2*ip.rhoCoeffs_;
192 }
193}
194
195
196template<class Specie, int PolySize>
198{
199 Specie::operator*=(s);
200}
201
202
203// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
204
205template<class Specie, int PolySize>
207(
210)
211{
212 Specie sp
213 (
214 static_cast<const Specie&>(ip1)
215 + static_cast<const Specie&>(ip2)
216 );
217
218 if (mag(sp.Y()) < SMALL)
219 {
221 (
222 sp,
223 ip1.rhoCoeffs_
224 );
225 }
226 else
227 {
228 const scalar Y1 = ip1.Y()/sp.Y();
229 const scalar Y2 = ip2.Y()/sp.Y();
230
231 return icoPolynomial<Specie, PolySize>
232 (
233 sp,
234 Y1*ip1.rhoCoeffs_ + Y2*ip2.rhoCoeffs_
235 );
236 }
237}
238
239
240template<class Specie, int PolySize>
242(
243 const scalar s,
244 const icoPolynomial<Specie, PolySize>& ip
245)
246{
247 return icoPolynomial<Specie, PolySize>
248 (
249 s*static_cast<const Specie&>(ip),
250 ip.rhoCoeffs_
251 );
252}
253
254
255template<class Specie, int PolySize>
257(
258 const icoPolynomial<Specie, PolySize>& ip1,
259 const icoPolynomial<Specie, PolySize>& ip2
260)
261{
262 Specie sp
263 (
264 static_cast<const Specie&>(ip1)
265 == static_cast<const Specie&>(ip2)
266 );
267
268 const scalar Y1 = ip1.Y()/sp.Y();
269 const scalar Y2 = ip2.Y()/sp.Y();
270
271 return icoPolynomial<Specie, PolySize>
272 (
273 sp,
274 Y2*ip2.rhoCoeffs_ - Y1*ip1.rhoCoeffs_
275 );
276}
277
278
279// ************************************************************************* //
Polynomial templated on size (order):
Definition: Polynomial.H:78
compactSpatialTensor S
The joint motion sub-space (3-DoF)
Definition: joint.H:132
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
Incompressible, polynomial form of equation of state, using a polynomial function for density.
autoPtr< icoPolynomial > clone() const
Construct and return a clone.
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
void operator*=(const scalar)
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & psi
const volScalarField & T
const volScalarField & Cv
Definition: EEqn.H:8
const volScalarField & Cp
Definition: EEqn.H:7
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))
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict