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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "hPolynomialThermo.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class EquationOfState, int PolySize>
34 (
35  const EquationOfState& pt,
36  const scalar Hf,
37  const scalar Sf,
38  const Polynomial<PolySize>& CpCoeffs,
39  const typename Polynomial<PolySize>::intPolyType& hCoeffs,
40  const Polynomial<PolySize>& sCoeffs
41 )
42 :
43  EquationOfState(pt),
44  Hf_(Hf),
45  Sf_(Sf),
46  CpCoeffs_(CpCoeffs),
47  hCoeffs_(hCoeffs),
48  sCoeffs_(sCoeffs)
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 template<class EquationOfState, int PolySize>
56 (
57  const word& name,
58  const hPolynomialThermo& pt
59 )
60 :
61  EquationOfState(name, pt),
62  Hf_(pt.Hf_),
63  Sf_(pt.Sf_),
64  CpCoeffs_(pt.CpCoeffs_),
65  hCoeffs_(pt.hCoeffs_),
66  sCoeffs_(pt.sCoeffs_)
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 template<class EquationOfState, int PolySize>
74 (
75  const scalar T
76 ) const
77 {
78  return T;
79 }
80 
81 
82 template<class EquationOfState, int PolySize>
84 (
85  const scalar p, const scalar T
86 ) const
87 {
88  return CpCoeffs_.value(T) + EquationOfState::Cp(p, T);
89 }
90 
91 
92 template<class EquationOfState, int PolySize>
94 (
95  const scalar p, const scalar T
96 ) const
97 {
98  return hCoeffs_.value(T) + EquationOfState::H(p, T);
99 }
100 
101 
102 template<class EquationOfState, int PolySize>
104 (
105  const scalar p, const scalar T
106 ) const
107 {
108  return Ha(p, T) - Hc();
109 }
110 
111 
112 template<class EquationOfState, int PolySize>
114 const
115 {
116  return Hf_;
117 }
118 
119 
120 template<class EquationOfState, int PolySize>
122 (
123  const scalar p,
124  const scalar T
125 ) const
126 {
127  return sCoeffs_.value(T) + EquationOfState::S(p, T);
128 }
129 
130 
131 template<class EquationOfState, int PolySize>
133 (
134  const scalar p,
135  const scalar T
136 ) const
137 {
138  return
139  (
140  hCoeffs_.derivative(T)
141  - T*sCoeffs_.derivative(T)
142  - sCoeffs_.value(T)
143  );
144 }
145 
146 
147 template<class EquationOfState, int PolySize>
149 (
150  const scalar p,
151  const scalar T
152 ) const
153 {
154  return
155  (
156  CpCoeffs_.derivative(T)
157  );
158 }
159 
160 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
161 
162 template<class EquationOfState, int PolySize>
164 (
166 )
167 {
168  scalar Y1 = this->Y();
169 
170  EquationOfState::operator+=(pt);
171 
172  if (mag(this->Y()) > SMALL)
173  {
174  Y1 /= this->Y();
175  const scalar Y2 = pt.Y()/this->Y();
176 
177  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
178  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
179  CpCoeffs_ = Y1*CpCoeffs_ + Y2*pt.CpCoeffs_;
180  hCoeffs_ = Y1*hCoeffs_ + Y2*pt.hCoeffs_;
181  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
182  }
183 }
184 
185 
186 template<class EquationOfState, int PolySize>
188 (
189  const scalar s
190 )
191 {
192  EquationOfState::operator*=(s);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
197 
198 template<class EquationOfState, int PolySize>
200 (
203 )
204 {
205  EquationOfState eofs = pt1;
206  eofs += pt2;
207 
208  if (mag(eofs.Y()) < SMALL)
209  {
211  (
212  eofs,
213  pt1.Hf_,
214  pt1.Sf_,
215  pt1.CpCoeffs_,
216  pt1.hCoeffs_,
217  pt1.sCoeffs_
218  );
219  }
220  {
221  const scalar Y1 = pt1.Y()/eofs.Y();
222  const scalar Y2 = pt2.Y()/eofs.Y();
223 
224  return hPolynomialThermo<EquationOfState, PolySize>
225  (
226  eofs,
227  Y1*pt1.Hf_ + Y2*pt2.Hf_,
228  Y1*pt1.Sf_ + Y2*pt2.Sf_,
229  Y1*pt1.CpCoeffs_ + Y2*pt2.CpCoeffs_,
230  Y1*pt1.hCoeffs_ + Y2*pt2.hCoeffs_,
231  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
232  );
233  }
234 }
235 
236 
237 template<class EquationOfState, int PolySize>
239 (
240  const scalar s,
241  const hPolynomialThermo<EquationOfState, PolySize>& pt
242 )
243 {
244  return hPolynomialThermo<EquationOfState, PolySize>
245  (
246  s*static_cast<const EquationOfState&>(pt),
247  pt.Hf_,
248  pt.Sf_,
249  pt.CpCoeffs_,
250  pt.hCoeffs_,
251  pt.sCoeffs_
252  );
253 }
254 
255 
256 template<class EquationOfState, int PolySize>
258 (
259  const hPolynomialThermo<EquationOfState, PolySize>& pt1,
260  const hPolynomialThermo<EquationOfState, PolySize>& pt2
261 )
262 {
263  EquationOfState eofs
264  (
265  static_cast<const EquationOfState&>(pt1)
266  == static_cast<const EquationOfState&>(pt2)
267  );
268 
269  const scalar Y1 = pt1.Y()/eofs.Y();
270  const scalar Y2 = pt2.Y()/eofs.Y();
271 
272  return hPolynomialThermo<EquationOfState, PolySize>
273  (
274  eofs,
275  Y2*pt2.Hf_ - Y1*pt1.Hf_,
276  Y2*pt2.Sf_ - Y1*pt1.Sf_,
277  Y2*pt2.CpCoeffs_ - Y1*pt1.CpCoeffs_,
278  Y2*pt2.hCoeffs_ - Y1*pt1.hCoeffs_,
279  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
280  );
281 }
282 
283 
284 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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 the temperature to be in the range Tlow_ to Thigh_.
Definition: hPolynomialThermoI.H:74
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::hPolynomialThermo::Hs
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hPolynomialThermoI.H:104
Foam::hPolynomialThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: hPolynomialThermoI.H:149
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::hPolynomialThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: hPolynomialThermoI.H:122
Foam::hPolynomialThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: hPolynomialThermoI.H:133
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:94
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:84
Ha
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:32
Foam::hPolynomialThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: hPolynomialThermoI.H:113