hPolynomialThermoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 template<class EquationOfState, int PolySize>
131 (
132  const scalar p,
133  const scalar T
134 ) const
135 {
136  return
137  (
138  hCoeffs_.derivative(T)
139  - T*sCoeffs_.derivative(T)
140  - sCoeffs_.value(T)
141  );
142 }
143 
144 
145 template<class EquationOfState, int PolySize>
147 (
148  const scalar p,
149  const scalar T
150 ) const
151 {
152  return
153  (
154  CpCoeffs_.derivative(T)
155  );
156 }
157 
158 
159 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
160 
161 template<class EquationOfState, int PolySize>
162 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator=
163 (
165 )
166 {
167  EquationOfState::operator=(pt);
168 
169  Hf_ = pt.Hf_;
170  Sf_ = pt.Sf_;
171  CpCoeffs_ = pt.CpCoeffs_;
172  hCoeffs_ = pt.hCoeffs_;
173  sCoeffs_ = pt.sCoeffs_;
174 }
175 
176 
177 template<class EquationOfState, int PolySize>
178 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator+=
179 (
181 )
182 {
183  scalar Y1 = this->Y();
184 
185  EquationOfState::operator+=(pt);
186 
187  if (mag(this->Y()) > small)
188  {
189  Y1 /= this->Y();
190  const scalar Y2 = pt.Y()/this->Y();
191 
192  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
193  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
194  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
195  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
196  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
197  }
198 }
199 
200 
201 template<class EquationOfState, int PolySize>
202 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
203 (
204  const scalar s
205 )
206 {
207  EquationOfState::operator*=(s);
208 }
209 
210 
211 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
212 
213 template<class EquationOfState, int PolySize>
215 (
218 )
219 {
220  EquationOfState eofs = pt1;
221  eofs += pt2;
222 
223  if (mag(eofs.Y()) < small)
224  {
226  (
227  eofs,
228  pt1.Hf_,
229  pt1.Sf_,
230  pt1.CpCoeffs_,
231  pt1.hCoeffs_,
232  pt1.sCoeffs_
233  );
234  }
235  {
236  const scalar Y1 = pt1.Y()/eofs.Y();
237  const scalar Y2 = pt2.Y()/eofs.Y();
238 
240  (
241  eofs,
242  Y1*pt1.Hf_ + Y2*pt2.Hf_,
243  Y1*pt1.Sf_ + Y2*pt2.Sf_,
244  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
245  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
246  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
247  );
248  }
249 }
250 
251 
252 template<class EquationOfState, int PolySize>
254 (
255  const scalar s,
257 )
258 {
260  (
261  s*static_cast<const EquationOfState&>(pt),
262  pt.Hf_,
263  pt.Sf_,
264  pt.CpCoeffs_,
265  pt.hCoeffs_,
266  pt.sCoeffs_
267  );
268 }
269 
270 
271 template<class EquationOfState, int PolySize>
273 (
276 )
277 {
278  EquationOfState eofs
279  (
280  static_cast<const EquationOfState&>(pt1)
281  == static_cast<const EquationOfState&>(pt2)
282  );
283 
284  const scalar Y1 = pt1.Y()/eofs.Y();
285  const scalar Y2 = pt2.Y()/eofs.Y();
286 
288  (
289  eofs,
290  Y2*pt2.Hf_ - Y1*pt1.Hf_,
291  Y2*pt2.Sf_ - Y1*pt1.Sf_,
292  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
293  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
294  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
295  );
296 }
297 
298 
299 // ************************************************************************* //
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
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].
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 dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].