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-2016 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) + EquationOfState::cp(p, 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) + EquationOfState::h(p, 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>
147 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator=
148 (
150 )
151 {
152  EquationOfState::operator=(pt);
153 
154  Hf_ = pt.Hf_;
155  Sf_ = pt.Sf_;
156  CpCoeffs_ = pt.CpCoeffs_;
157  hCoeffs_ = pt.hCoeffs_;
158  sCoeffs_ = pt.sCoeffs_;
159 }
160 
161 
162 template<class EquationOfState, int PolySize>
163 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator+=
164 (
166 )
167 {
168  scalar molr1 = this->nMoles();
169 
170  EquationOfState::operator+=(pt);
171 
172  molr1 /= this->nMoles();
173  scalar molr2 = pt.nMoles()/this->nMoles();
174 
175  Hf_ = molr1*Hf_ + molr2*pt.Hf_;
176  Sf_ = molr1*Sf_ + molr2*pt.Sf_;
177  CpCoeffs_ = molr1*CpCoeffs_ + molr2*pt.CpCoeffs_;
178  hCoeffs_ = molr1*hCoeffs_ + molr2*pt.hCoeffs_;
179  sCoeffs_ = molr1*sCoeffs_ + molr2*pt.sCoeffs_;
180 }
181 
182 
183 template<class EquationOfState, int PolySize>
184 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator-=
185 (
187 )
188 {
189  scalar molr1 = this->nMoles();
190 
191  EquationOfState::operator-=(pt);
192 
193  molr1 /= this->nMoles();
194  scalar molr2 = pt.nMoles()/this->nMoles();
195 
196  Hf_ = molr1*Hf_ - molr2*pt.Hf_;
197  Sf_ = molr1*Sf_ - molr2*pt.Sf_;
198  CpCoeffs_ = molr1*CpCoeffs_ - molr2*pt.CpCoeffs_;
199  hCoeffs_ = molr1*hCoeffs_ - molr2*pt.hCoeffs_;
200  sCoeffs_ = molr1*sCoeffs_ - molr2*pt.sCoeffs_;
201 }
202 
203 
204 template<class EquationOfState, int PolySize>
205 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
206 (
207  const scalar s
208 )
209 {
210  EquationOfState::operator*=(s);
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
215 
216 template<class EquationOfState, int PolySize>
218 (
221 )
222 {
223  EquationOfState eofs = pt1;
224  eofs += pt2;
225 
226  scalar molr1 = pt1.nMoles()/eofs.nMoles();
227  scalar molr2 = pt2.nMoles()/eofs.nMoles();
228 
230  (
231  eofs,
232  molr1*pt1.Hf_ + molr2*pt2.Hf_,
233  molr1*pt1.Sf_ + molr2*pt2.Sf_,
234  molr1*pt1.CpCoeffs_ + molr2*pt2.CpCoeffs_,
235  molr1*pt1.hCoeffs_ + molr2*pt2.hCoeffs_,
236  molr1*pt1.sCoeffs_ + molr2*pt2.sCoeffs_
237  );
238 }
239 
240 
241 template<class EquationOfState, int PolySize>
243 (
246 )
247 {
248  EquationOfState eofs = pt1;
249  eofs -= pt2;
250 
251  scalar molr1 = pt1.nMoles()/eofs.nMoles();
252  scalar molr2 = pt2.nMoles()/eofs.nMoles();
253 
255  (
256  eofs,
257  molr1*pt1.Hf_ - molr2*pt2.Hf_,
258  molr1*pt1.Sf_ - molr2*pt2.Sf_,
259  molr1*pt1.CpCoeffs_ - molr2*pt2.CpCoeffs_,
260  molr1*pt1.hCoeffs_ - molr2*pt2.hCoeffs_,
261  molr1*pt1.sCoeffs_ - molr2*pt2.sCoeffs_
262  );
263 }
264 
265 
266 template<class EquationOfState, int PolySize>
268 (
269  const scalar s,
271 )
272 {
274  (
275  s*static_cast<const EquationOfState&>(pt),
276  pt.Hf_,
277  pt.Sf_,
278  pt.CpCoeffs_,
279  pt.hCoeffs_,
280  pt.sCoeffs_
281  );
282 }
283 
284 
285 template<class EquationOfState, int PolySize>
287 (
290 )
291 {
292  return pt2 - pt1;
293 }
294 
295 
296 // ************************************************************************* //
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
scalar hc() const
Chemical enthalpy [J/kmol].
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Thermodynamics package templated on the equation of state, using polynomial functions for cp...
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
const volScalarField & cp
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kmol].