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-2019 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  scalar Y1 = this->Y();
168 
169  EquationOfState::operator+=(pt);
170 
171  if (mag(this->Y()) > small)
172  {
173  Y1 /= this->Y();
174  const scalar Y2 = pt.Y()/this->Y();
175 
176  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
177  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
178  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
179  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
180  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
181  }
182 }
183 
184 
185 template<class EquationOfState, int PolySize>
186 inline void Foam::hPolynomialThermo<EquationOfState, PolySize>::operator*=
187 (
188  const scalar s
189 )
190 {
191  EquationOfState::operator*=(s);
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
196 
197 template<class EquationOfState, int PolySize>
199 (
202 )
203 {
204  EquationOfState eofs = pt1;
205  eofs += pt2;
206 
207  if (mag(eofs.Y()) < small)
208  {
210  (
211  eofs,
212  pt1.Hf_,
213  pt1.Sf_,
214  pt1.CpCoeffs_,
215  pt1.hCoeffs_,
216  pt1.sCoeffs_
217  );
218  }
219  {
220  const scalar Y1 = pt1.Y()/eofs.Y();
221  const scalar Y2 = pt2.Y()/eofs.Y();
222 
224  (
225  eofs,
226  Y1*pt1.Hf_ + Y2*pt2.Hf_,
227  Y1*pt1.Sf_ + Y2*pt2.Sf_,
228  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
229  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
230  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
231  );
232  }
233 }
234 
235 
236 template<class EquationOfState, int PolySize>
238 (
239  const scalar s,
241 )
242 {
244  (
245  s*static_cast<const EquationOfState&>(pt),
246  pt.Hf_,
247  pt.Sf_,
248  pt.CpCoeffs_,
249  pt.hCoeffs_,
250  pt.sCoeffs_
251  );
252 }
253 
254 
255 template<class EquationOfState, int PolySize>
257 (
260 )
261 {
262  EquationOfState eofs
263  (
264  static_cast<const EquationOfState&>(pt1)
265  == static_cast<const EquationOfState&>(pt2)
266  );
267 
268  const scalar Y1 = pt1.Y()/eofs.Y();
269  const scalar Y2 = pt2.Y()/eofs.Y();
270 
272  (
273  eofs,
274  Y2*pt2.Hf_ - Y1*pt1.Hf_,
275  Y2*pt2.Sf_ - Y1*pt1.Sf_,
276  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
277  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
278  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
279  );
280 }
281 
282 
283 // ************************************************************************* //
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
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
mesh Sf()
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].