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-2020 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,
84  const scalar T
85 ) const
86 {
87  return CpCoeffs_.value(T) + EquationOfState::Cp(p, T);
88 }
89 
90 
91 template<class EquationOfState, int PolySize>
93 (
94  const scalar p,
95  const scalar T
96 ) const
97 {
98  return Ha(p, T) - Hf();
99 }
100 
101 
102 template<class EquationOfState, int PolySize>
104 (
105  const scalar p,
106  const scalar T
107 ) const
108 {
109  return hCoeffs_.value(T) + EquationOfState::H(p, T);
110 }
111 
112 
113 template<class EquationOfState, int PolySize>
115 const
116 {
117  return Hf_;
118 }
119 
120 
121 template<class EquationOfState, int PolySize>
123 (
124  const scalar p,
125  const scalar T
126 ) const
127 {
128  return sCoeffs_.value(T) + EquationOfState::Sp(p, T);
129 }
130 
131 
132 template<class EquationOfState, int PolySize>
134 (
135  const scalar T
136 ) const
137 {
138  return hCoeffs_.value(T) - sCoeffs_.value(T)*T;
139 }
140 
141 
142 template<class EquationOfState, int PolySize>
144 (
145  const scalar p,
146  const scalar T
147 ) const
148 {
149  return
150  (
151  CpCoeffs_.derivative(T)
152  );
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157 
158 template<class EquationOfState, int PolySize>
159 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator+=
160 (
162 )
163 {
164  scalar Y1 = this->Y();
165 
166  EquationOfState::operator+=(pt);
167 
168  if (mag(this->Y()) > small)
169  {
170  Y1 /= this->Y();
171  const scalar Y2 = pt.Y()/this->Y();
172 
173  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
174  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
175  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
176  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
177  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
178  }
179 }
180 
181 
182 template<class EquationOfState, int PolySize>
183 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
184 (
185  const scalar s
186 )
187 {
188  EquationOfState::operator*=(s);
189 }
190 
191 
192 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
193 
194 template<class EquationOfState, int PolySize>
196 (
199 )
200 {
201  EquationOfState eofs = pt1;
202  eofs += pt2;
203 
204  if (mag(eofs.Y()) < small)
205  {
207  (
208  eofs,
209  pt1.Hf_,
210  pt1.Sf_,
211  pt1.CpCoeffs_,
212  pt1.hCoeffs_,
213  pt1.sCoeffs_
214  );
215  }
216  {
217  const scalar Y1 = pt1.Y()/eofs.Y();
218  const scalar Y2 = pt2.Y()/eofs.Y();
219 
221  (
222  eofs,
223  Y1*pt1.Hf_ + Y2*pt2.Hf_,
224  Y1*pt1.Sf_ + Y2*pt2.Sf_,
225  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
226  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
227  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
228  );
229  }
230 }
231 
232 
233 template<class EquationOfState, int PolySize>
235 (
236  const scalar s,
238 )
239 {
241  (
242  s*static_cast<const EquationOfState&>(pt),
243  pt.Hf_,
244  pt.Sf_,
245  pt.CpCoeffs_,
246  pt.hCoeffs_,
247  pt.sCoeffs_
248  );
249 }
250 
251 
252 template<class EquationOfState, int PolySize>
254 (
257 )
258 {
259  EquationOfState eofs
260  (
261  static_cast<const EquationOfState&>(pt1)
262  == static_cast<const EquationOfState&>(pt2)
263  );
264 
265  const scalar Y1 = pt1.Y()/eofs.Y();
266  const scalar Y2 = pt2.Y()/eofs.Y();
267 
269  (
270  eofs,
271  Y2*pt2.Hf_ - Y1*pt1.Hf_,
272  Y2*pt2.Sf_ - Y1*pt1.Sf_,
273  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
274  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
275  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
276  );
277 }
278 
279 
280 // ************************************************************************* //
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].
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
Enthalpy based thermodynamics package using a polynomial function of temperature for the constant hea...
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 Hf() const
Enthalpy of formation [J/kg].
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].