ePolynomialThermoI.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) 2020-2023 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 "ePolynomialThermo.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class EquationOfState, int PolySize>
32 (
33  const EquationOfState& pt,
34  const scalar Hf,
35  const scalar Sf,
36  const Polynomial<PolySize>& CvCoeffs,
37  const typename Polynomial<PolySize>::intPolyType& eCoeffs,
38  const Polynomial<PolySize>& sCoeffs
39 )
40 :
41  EquationOfState(pt),
42  Hf_(Hf),
43  Sf_(Sf),
44  CvCoeffs_(CvCoeffs),
45  eCoeffs_(eCoeffs),
46  sCoeffs_(sCoeffs)
47 {}
48 
49 
50 template<class EquationOfState, int PolySize>
52 (
53  const word& name,
54  const ePolynomialThermo& pt
55 )
56 :
57  EquationOfState(name, pt),
58  Hf_(pt.Hf_),
59  Sf_(pt.Sf_),
60  CvCoeffs_(pt.CvCoeffs_),
61  eCoeffs_(pt.eCoeffs_),
62  sCoeffs_(pt.sCoeffs_)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
68 template<class EquationOfState, int PolySize>
70 (
71  const scalar T
72 ) const
73 {
74  return T;
75 }
76 
77 
78 template<class EquationOfState, int PolySize>
80 (
81  const scalar p,
82  const scalar T
83 ) const
84 {
85  return CvCoeffs_.value(T) + EquationOfState::Cv(p, T);
86 }
87 
88 
89 template<class EquationOfState, int PolySize>
91 (
92  const scalar p,
93  const scalar T
94 ) const
95 {
96  return eCoeffs_.value(T) + EquationOfState::E(p, T);
97 }
98 
99 
100 template<class EquationOfState, int PolySize>
102 (
103  const scalar p,
104  const scalar T
105 ) const
106 {
107  return Es(p, T) + Hf();
108 }
109 
110 
111 template<class EquationOfState, int PolySize>
113 const
114 {
115  return Hf_;
116 }
117 
118 
119 template<class EquationOfState, int PolySize>
121 (
122  const scalar p,
123  const scalar T
124 ) const
125 {
126  return sCoeffs_.value(T) + EquationOfState::Sv(p, T);
127 }
128 
129 
130 template<class EquationOfState, int PolySize>
132 (
133  const scalar T
134 ) const
135 {
136  return
137  eCoeffs_.value(T) + Hf() + Pstd/EquationOfState::rho(Pstd, T)
138  - S(Pstd, T)*T;
139 }
140 
141 
142 template<class EquationOfState, int PolySize>
144 (
145  const scalar p,
146  const scalar T
147 ) const
148 {
150  return CvCoeffs_.derivative(T); // + EquationOfState::dCpdT
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
155 
156 template<class EquationOfState, int PolySize>
158 (
160 )
161 {
162  scalar Y1 = this->Y();
163 
164  EquationOfState::operator+=(pt);
165 
166  if (mag(this->Y()) > small)
167  {
168  Y1 /= this->Y();
169  const scalar Y2 = pt.Y()/this->Y();
170 
171  Hf_ = Y1*Hf_ + Y2*pt.Hf_;
172  Sf_ = Y1*Sf_ + Y2*pt.Sf_;
173  CvCoeffs_ = Y1*CvCoeffs_ + Y2*pt.CvCoeffs_;
174  eCoeffs_ = Y1*eCoeffs_ + Y2*pt.eCoeffs_;
175  sCoeffs_ = Y1*sCoeffs_ + Y2*pt.sCoeffs_;
176  }
177 }
178 
179 
180 template<class EquationOfState, int PolySize>
182 (
183  const scalar s
184 )
185 {
186  EquationOfState::operator*=(s);
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
191 
192 template<class EquationOfState, int PolySize>
194 (
197 )
198 {
199  EquationOfState eofs = pt1;
200  eofs += pt2;
201 
202  if (mag(eofs.Y()) < small)
203  {
205  (
206  eofs,
207  pt1.Hf_,
208  pt1.Sf_,
209  pt1.CvCoeffs_,
210  pt1.eCoeffs_,
211  pt1.sCoeffs_
212  );
213  }
214  {
215  const scalar Y1 = pt1.Y()/eofs.Y();
216  const scalar Y2 = pt2.Y()/eofs.Y();
217 
218  return ePolynomialThermo<EquationOfState, PolySize>
219  (
220  eofs,
221  Y1*pt1.Hf_ + Y2*pt2.Hf_,
222  Y1*pt1.Sf_ + Y2*pt2.Sf_,
223  Y1*pt1.CvCoeffs_ + Y2*pt2.CvCoeffs_,
224  Y1*pt1.eCoeffs_ + Y2*pt2.eCoeffs_,
225  Y1*pt1.sCoeffs_ + Y2*pt2.sCoeffs_
226  );
227  }
228 }
229 
230 
231 template<class EquationOfState, int PolySize>
233 (
234  const scalar s,
235  const ePolynomialThermo<EquationOfState, PolySize>& pt
236 )
237 {
238  return ePolynomialThermo<EquationOfState, PolySize>
239  (
240  s*static_cast<const EquationOfState&>(pt),
241  pt.Hf_,
242  pt.Sf_,
243  pt.CvCoeffs_,
244  pt.eCoeffs_,
245  pt.sCoeffs_
246  );
247 }
248 
249 
250 template<class EquationOfState, int PolySize>
252 (
253  const ePolynomialThermo<EquationOfState, PolySize>& pt1,
254  const ePolynomialThermo<EquationOfState, PolySize>& pt2
255 )
256 {
257  EquationOfState eofs
258  (
259  static_cast<const EquationOfState&>(pt1)
260  == static_cast<const EquationOfState&>(pt2)
261  );
262 
263  const scalar Y1 = pt1.Y()/eofs.Y();
264  const scalar Y2 = pt2.Y()/eofs.Y();
265 
266  return ePolynomialThermo<EquationOfState, PolySize>
267  (
268  eofs,
269  Y2*pt2.Hf_ - Y1*pt1.Hf_,
270  Y2*pt2.Sf_ - Y1*pt1.Sf_,
271  Y2*pt2.CvCoeffs_ - Y1*pt1.CvCoeffs_,
272  Y2*pt2.eCoeffs_ - Y1*pt1.eCoeffs_,
273  Y2*pt2.sCoeffs_ - Y1*pt1.sCoeffs_
274  );
275 }
276 
277 
278 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
scalar Es(const scalar p, const scalar T) const
Definition: HtoEthermo.H:11
Internal energy based thermodynamics package using a polynomial function of temperature for the const...
scalar Hf() const
Enthalpy of formation [J/kg].
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/kg/K].
scalar limit(const scalar) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar Pstd
Standard pressure.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p
PtrList< volScalarField > & Y