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-2017 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 word& name,
56  const hPolynomialThermo& pt
57 )
58 :
59  EquationOfState(name, pt),
60  Hf_(pt.Hf_),
61  Sf_(pt.Sf_),
62  CpCoeffs_(pt.CpCoeffs_),
63  hCoeffs_(pt.hCoeffs_),
64  sCoeffs_(pt.sCoeffs_)
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class EquationOfState, int PolySize>
72 (
73  const scalar T
74 ) const
75 {
76  return T;
77 }
78 
79 
80 template<class EquationOfState, int PolySize>
82 (
83  const scalar p, const scalar T
84 ) const
85 {
86  return CpCoeffs_.value(T) + EquationOfState::Cp(p, T);
87 }
88 
89 
90 template<class EquationOfState, int PolySize>
92 (
93  const scalar p, const scalar T
94 ) const
95 {
96  return hCoeffs_.value(T) + EquationOfState::H(p, T);
97 }
98 
99 
100 template<class EquationOfState, int PolySize>
102 (
103  const scalar p, const scalar T
104 ) const
105 {
106  return Ha(p, T) - Hc();
107 }
108 
109 
110 template<class EquationOfState, int PolySize>
112 const
113 {
114  return Hf_;
115 }
116 
117 
118 template<class EquationOfState, int PolySize>
120 (
121  const scalar p,
122  const scalar T
123 ) const
124 {
125  return sCoeffs_.value(T) + EquationOfState::S(p, T);
126 }
127 
128 
129 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
130 
131 template<class EquationOfState, int PolySize>
132 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator=
133 (
135 )
136 {
137  EquationOfState::operator=(pt);
138 
139  Hf_ = pt.Hf_;
140  Sf_ = pt.Sf_;
141  CpCoeffs_ = pt.CpCoeffs_;
142  hCoeffs_ = pt.hCoeffs_;
143  sCoeffs_ = pt.sCoeffs_;
144 }
145 
146 
147 template<class EquationOfState, int PolySize>
148 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator+=
149 (
151 )
152 {
153  scalar Y1 = this->Y();
154 
155  EquationOfState::operator+=(pt);
156 
157  if (mag(this->Y()) > SMALL)
158  {
159  Y1 /= this->Y();
160  const scalar Y2 = pt.Y()/this->Y();
161 
162  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
163  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
164  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
165  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
166  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
167  }
168 }
169 
170 
171 template<class EquationOfState, int PolySize>
172 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
173 (
174  const scalar s
175 )
176 {
177  EquationOfState::operator*=(s);
178 }
179 
180 
181 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
182 
183 template<class EquationOfState, int PolySize>
185 (
188 )
189 {
190  EquationOfState eofs = pt1;
191  eofs += pt2;
192 
193  if (mag(eofs.Y()) < SMALL)
194  {
196  (
197  eofs,
198  pt1.Hf_,
199  pt1.Sf_,
200  pt1.CpCoeffs_,
201  pt1.hCoeffs_,
202  pt1.sCoeffs_
203  );
204  }
205  {
206  const scalar Y1 = pt1.Y()/eofs.Y();
207  const scalar Y2 = pt2.Y()/eofs.Y();
208 
210  (
211  eofs,
212  Y1*pt1.Hf_ + Y2*pt2.Hf_,
213  Y1*pt1.Sf_ + Y2*pt2.Sf_,
214  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
215  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
216  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
217  );
218  }
219 }
220 
221 
222 template<class EquationOfState, int PolySize>
224 (
225  const scalar s,
227 )
228 {
230  (
231  s*static_cast<const EquationOfState&>(pt),
232  pt.Hf_,
233  pt.Sf_,
234  pt.CpCoeffs_,
235  pt.hCoeffs_,
236  pt.sCoeffs_
237  );
238 }
239 
240 
241 template<class EquationOfState, int PolySize>
243 (
246 )
247 {
248  EquationOfState eofs
249  (
250  static_cast<const EquationOfState&>(pt1)
251  == static_cast<const EquationOfState&>(pt2)
252  );
253 
254  const scalar Y1 = pt1.Y()/eofs.Y();
255  const scalar Y2 = pt2.Y()/eofs.Y();
256 
258  (
259  eofs,
260  Y2*pt2.Hf_ - Y1*pt1.Hf_,
261  Y2*pt2.Sf_ - Y1*pt1.Sf_,
262  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
263  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
264  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
265  );
266 }
267 
268 
269 // ************************************************************************* //
PtrList< volScalarField > & Y1
Definition: YEqns.H:8
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Thermodynamics package templated on the equation of state, using polynomial functions for cp...
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
PtrList< volScalarField > & Y2
Definition: YEqns.H:9
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar Hc() const
Chemical enthalpy [J/kg].
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].