hPolynomialThermoI.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 | Copyright (C) 2011-2015 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 "hPolynomialThermo.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class EquationOfState, int PolySize>
32 (
33  const EquationOfState& pt,
34  const scalar Hf,
35  const scalar Sf,
36  const Polynomial<PolySize>& CpCoeffs,
37  const typename Polynomial<PolySize>::intPolyType& hCoeffs,
38  const Polynomial<PolySize>& sCoeffs
39 )
40 :
41  EquationOfState(pt),
42  Hf_(Hf),
43  Sf_(Sf),
44  CpCoeffs_(CpCoeffs),
45  hCoeffs_(hCoeffs),
46  sCoeffs_(sCoeffs)
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class EquationOfState, int PolySize>
54 (
55  const hPolynomialThermo& pt
56 )
57 :
58  EquationOfState(pt),
59  Hf_(pt.Hf_),
60  Sf_(pt.Sf_),
61  CpCoeffs_(pt.CpCoeffs_),
62  hCoeffs_(pt.hCoeffs_),
63  sCoeffs_(pt.sCoeffs_)
64 {}
65 
66 
67 template<class EquationOfState, int PolySize>
69 (
70  const word& name,
71  const hPolynomialThermo& pt
72 )
73 :
74  EquationOfState(name, pt),
75  Hf_(pt.Hf_),
76  Sf_(pt.Sf_),
77  CpCoeffs_(pt.CpCoeffs_),
78  hCoeffs_(pt.hCoeffs_),
79  sCoeffs_(pt.sCoeffs_)
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class EquationOfState, int PolySize>
87 (
88  const scalar T
89 ) const
90 {
91  return T;
92 }
93 
94 
95 template<class EquationOfState, int PolySize>
97 (
98  const scalar p, const scalar T
99 ) const
100 {
101  return CpCoeffs_.value(T);
102 }
103 
104 
105 template<class EquationOfState, int PolySize>
107 (
108  const scalar p, const scalar T
109 ) const
110 {
111  return hCoeffs_.value(T);
112 }
113 
114 
115 template<class EquationOfState, int PolySize>
117 (
118  const scalar p, const scalar T
119 ) const
120 {
121  return ha(p, T) - hc();
122 }
123 
124 
125 template<class EquationOfState, int PolySize>
127 const
128 {
129  return Hf_;
130 }
131 
132 
133 template<class EquationOfState, int PolySize>
135 (
136  const scalar p,
137  const scalar T
138 ) const
139 {
140  return sCoeffs_.value(T) + EquationOfState::s(p, T);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
145 
146 template<class EquationOfState, int PolySize>
149 (
151 )
152 {
153  EquationOfState::operator=(pt);
154 
155  Hf_ = pt.Hf_;
156  Sf_ = pt.Sf_;
157  CpCoeffs_ = pt.CpCoeffs_;
158  hCoeffs_ = pt.hCoeffs_;
159  sCoeffs_ = pt.sCoeffs_;
160 
161  return *this;
162 }
163 
164 
165 template<class EquationOfState, int PolySize>
167 (
169 )
170 {
171  scalar molr1 = this->nMoles();
172 
173  EquationOfState::operator+=(pt);
174 
175  molr1 /= this->nMoles();
176  scalar molr2 = pt.nMoles()/this->nMoles();
177 
178  Hf_ = molr1*Hf_ + molr2*pt.Hf_;
179  Sf_ = molr1*Sf_ + molr2*pt.Sf_;
180  CpCoeffs_ = molr1*CpCoeffs_ + molr2*pt.CpCoeffs_;
181  hCoeffs_ = molr1*hCoeffs_ + molr2*pt.hCoeffs_;
182  sCoeffs_ = molr1*sCoeffs_ + molr2*pt.sCoeffs_;
183 }
184 
185 
186 template<class EquationOfState, int PolySize>
188 (
190 )
191 {
192  scalar molr1 = this->nMoles();
193 
194  EquationOfState::operator-=(pt);
195 
196  molr1 /= this->nMoles();
197  scalar molr2 = pt.nMoles()/this->nMoles();
198 
199  Hf_ = molr1*Hf_ - molr2*pt.Hf_;
200  Sf_ = molr1*Sf_ - molr2*pt.Sf_;
201  CpCoeffs_ = molr1*CpCoeffs_ - molr2*pt.CpCoeffs_;
202  hCoeffs_ = molr1*hCoeffs_ - molr2*pt.hCoeffs_;
203  sCoeffs_ = molr1*sCoeffs_ - molr2*pt.sCoeffs_;
204 }
205 
206 
207 template<class EquationOfState, int PolySize>
209 (
210  const scalar s
211 )
212 {
213  EquationOfState::operator*=(s);
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
218 
219 template<class EquationOfState, int PolySize>
221 (
224 )
225 {
226  EquationOfState eofs = pt1;
227  eofs += pt2;
228 
229  scalar molr1 = pt1.nMoles()/eofs.nMoles();
230  scalar molr2 = pt2.nMoles()/eofs.nMoles();
231 
233  (
234  eofs,
235  molr1*pt1.Hf_ + molr2*pt2.Hf_,
236  molr1*pt1.Sf_ + molr2*pt2.Sf_,
237  molr1*pt1.CpCoeffs_ + molr2*pt2.CpCoeffs_,
238  molr1*pt1.hCoeffs_ + molr2*pt2.hCoeffs_,
239  molr1*pt1.sCoeffs_ + molr2*pt2.sCoeffs_
240  );
241 }
242 
243 
244 template<class EquationOfState, int PolySize>
246 (
247  const hPolynomialThermo<EquationOfState, PolySize>& pt1,
248  const hPolynomialThermo<EquationOfState, PolySize>& pt2
249 )
250 {
251  EquationOfState eofs = pt1;
252  eofs -= pt2;
253 
254  scalar molr1 = pt1.nMoles()/eofs.nMoles();
255  scalar molr2 = pt2.nMoles()/eofs.nMoles();
256 
257  return hPolynomialThermo<EquationOfState, PolySize>
258  (
259  eofs,
260  molr1*pt1.Hf_ - molr2*pt2.Hf_,
261  molr1*pt1.Sf_ - molr2*pt2.Sf_,
262  molr1*pt1.CpCoeffs_ - molr2*pt2.CpCoeffs_,
263  molr1*pt1.hCoeffs_ - molr2*pt2.hCoeffs_,
264  molr1*pt1.sCoeffs_ - molr2*pt2.sCoeffs_
265  );
266 }
267 
268 
269 template<class EquationOfState, int PolySize>
271 (
272  const scalar s,
273  const hPolynomialThermo<EquationOfState, PolySize>& pt
274 )
275 {
276  return hPolynomialThermo<EquationOfState, PolySize>
277  (
278  s*static_cast<const EquationOfState&>(pt),
279  pt.Hf_,
280  pt.Sf_,
281  pt.CpCoeffs_,
282  pt.hCoeffs_,
283  pt.sCoeffs_
284  );
285 }
286 
287 
288 template<class EquationOfState, int PolySize>
290 (
291  const hPolynomialThermo<EquationOfState, PolySize>& pt1,
292  const hPolynomialThermo<EquationOfState, PolySize>& pt2
293 )
294 {
295  return pt2 - pt1;
296 }
297 
298 
299 // ************************************************************************* //
Foam::hPolynomialThermo::s
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: hPolynomialThermoI.H:135
p
p
Definition: pEqn.H:62
Foam::hPolynomialThermo::limit
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hPolynomialThermoI.H:87
Foam::hPolynomialThermo::hc
scalar hc() const
Chemical enthalpy [J/kmol].
Definition: hPolynomialThermoI.H:126
Foam::hPolynomialThermo::ha
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
Definition: hPolynomialThermoI.H:107
Foam::hPolynomialThermo::hs
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].
Definition: hPolynomialThermoI.H:117
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
s
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))
Foam::hPolynomialThermo
Thermodynamics package templated on the equation of state, using polynomial functions for cp,...
Definition: hPolynomialThermo.H:52
Foam::Polynomial< PolySize >
hPolynomialThermo.H
Foam::hPolynomialThermo::hPolynomialThermo
hPolynomialThermo(const EquationOfState &pt, const scalar Hf, const scalar Sf, const Polynomial< PolySize > &CpCoeffs, const typename Polynomial< PolySize >::intPolyType &hCoeffs, const Polynomial< PolySize > &sCoeffs)
Construct from components.
Definition: hPolynomialThermoI.H:32
Foam::hPolynomialThermo::cp
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
Definition: hPolynomialThermoI.H:97
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47