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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "hPolynomialThermo.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class EquationOfState, int PolySize>
35 (
36  const EquationOfState& pt,
37  const scalar Hf,
38  const scalar Sf,
39  const Polynomial<PolySize>& CpCoeffs,
40  const typename Polynomial<PolySize>::intPolyType& hCoeffs,
41  const Polynomial<PolySize>& sCoeffs
42 )
43 :
44  EquationOfState(pt),
45  Hf_(Hf),
46  Sf_(Sf),
47  CpCoeffs_(CpCoeffs),
48  hCoeffs_(hCoeffs),
49  sCoeffs_(sCoeffs)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 template<class EquationOfState, int PolySize>
57 (
58  const word& name,
59  const hPolynomialThermo& pt
60 )
61 :
62  EquationOfState(name, pt),
63  Hf_(pt.Hf_),
64  Sf_(pt.Sf_),
65  CpCoeffs_(pt.CpCoeffs_),
66  hCoeffs_(pt.hCoeffs_),
67  sCoeffs_(pt.sCoeffs_)
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class EquationOfState, int PolySize>
75 (
76  const scalar T
77 ) const
78 {
79  return T;
80 }
81 
82 
83 template<class EquationOfState, int PolySize>
85 (
86  const scalar p, const scalar T
87 ) const
88 {
89  return CpCoeffs_.value(T) + EquationOfState::Cp(p, T);
90 }
91 
92 
93 template<class EquationOfState, int PolySize>
95 (
96  const scalar p, const scalar T
97 ) const
98 {
99  return hCoeffs_.value(T) + EquationOfState::H(p, T);
100 }
101 
102 
103 template<class EquationOfState, int PolySize>
105 (
106  const scalar p, const scalar T
107 ) const
108 {
109  return Ha(p, T) - Hc();
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::S(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 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
156 
157 template<class EquationOfState, int PolySize>
159 (
161 )
162 {
163  scalar Y1 = this->Y();
164 
165  EquationOfState::operator+=(pt);
166 
167  if (mag(this->Y()) > SMALL)
168  {
169  Y1 /= this->Y();
170  const scalar Y2 = pt.Y()/this->Y();
171 
172  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
173  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
174  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
175  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
176  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
177  }
178 }
179 
180 
181 template<class EquationOfState, int PolySize>
183 (
184  const scalar s
185 )
186 {
187  EquationOfState::operator*=(s);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
192 
193 template<class EquationOfState, int PolySize>
195 (
198 )
199 {
200  EquationOfState eofs = pt1;
201  eofs += pt2;
202 
203  if (mag(eofs.Y()) < SMALL)
204  {
206  (
207  eofs,
208  pt1.Hf_,
209  pt1.Sf_,
210  pt1.CpCoeffs_,
211  pt1.hCoeffs_,
212  pt1.sCoeffs_
213  );
214  }
215  {
216  const scalar Y1 = pt1.Y()/eofs.Y();
217  const scalar Y2 = pt2.Y()/eofs.Y();
218 
219  return hPolynomialThermo<EquationOfState, PolySize>
220  (
221  eofs,
222  Y1*pt1.Hf_ + Y2*pt2.Hf_,
223  Y1*pt1.Sf_ + Y2*pt2.Sf_,
224  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
225  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
226  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
227  );
228  }
229 }
230 
231 
232 template<class EquationOfState, int PolySize>
234 (
235  const scalar s,
236  const hPolynomialThermo<EquationOfState, PolySize>& pt
237 )
238 {
239  return hPolynomialThermo<EquationOfState, PolySize>
240  (
241  s*static_cast<const EquationOfState&>(pt),
242  pt.Hf_,
243  pt.Sf_,
244  pt.CpCoeffs_,
245  pt.hCoeffs_,
246  pt.sCoeffs_
247  );
248 }
249 
250 
251 template<class EquationOfState, int PolySize>
253 (
254  const hPolynomialThermo<EquationOfState, PolySize>& pt1,
255  const hPolynomialThermo<EquationOfState, PolySize>& pt2
256 )
257 {
258  EquationOfState eofs
259  (
260  static_cast<const EquationOfState&>(pt1)
261  == static_cast<const EquationOfState&>(pt2)
262  );
263 
264  const scalar Y1 = pt1.Y()/eofs.Y();
265  const scalar Y2 = pt2.Y()/eofs.Y();
266 
267  return hPolynomialThermo<EquationOfState, PolySize>
268  (
269  eofs,
270  Y2*pt2.Hf_ - Y1*pt1.Hf_,
271  Y2*pt2.Sf_ - Y1*pt1.Sf_,
272  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
273  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
274  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
275  );
276 }
277 
278 
279 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::hPolynomialThermo::limit
scalar limit(const scalar) const
Limit temperature to be within the range.
Definition: hPolynomialThermoI.H:75
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::hPolynomialThermo::Gstd
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
Definition: hPolynomialThermoI.H:134
Foam::hPolynomialThermo::Hs
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hPolynomialThermoI.H:105
Foam::hPolynomialThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: hPolynomialThermoI.H:144
Foam::hPolynomialThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: hPolynomialThermoI.H:123
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::hPolynomialThermo::Ha
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
Definition: hPolynomialThermoI.H:95
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::hPolynomialThermo
Thermodynamics package templated on the equation of state, using polynomial functions for cp,...
Definition: hPolynomialThermo.H:104
hPolynomialThermo.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::hPolynomialThermo::Cp
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: hPolynomialThermoI.H:85
Ha
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:32
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::hPolynomialThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: hPolynomialThermoI.H:114