ePowerThermoI.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 "ePowerThermo.H"
27 #include "specie.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class EquationOfState>
33 (
34  const scalar T
35 ) const
36 {
37  if (T < 0)
38  {
40  << "attempt to evaluate ePowerThermo<EquationOfState>"
41  " for negative temperature " << T
42  << abort(FatalError);
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class EquationOfState>
51 (
52  const EquationOfState& st,
53  const scalar c0,
54  const scalar n0,
55  const scalar Tref,
56  const scalar Hf
57 )
58 :
59  EquationOfState(st),
60  c0_(c0),
61  n0_(n0),
62  Tref_(Tref),
63  Hf_(Hf)
64 {}
65 
66 
67 template<class EquationOfState>
69 (
70  const word& name,
71  const ePowerThermo& jt
72 )
73 :
74  EquationOfState(name, jt),
75  c0_(jt.c0_),
76  n0_(jt.n0_),
77  Tref_(jt.Tref_),
78  Hf_(jt.Hf_)
79 {}
80 
81 
82 template<class EquationOfState>
85 {
87  (
89  );
90 }
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 template<class EquationOfState>
97 (
98  const scalar T
99 ) const
100 {
101  return T;
102 }
103 
104 
105 template<class EquationOfState>
107 (
108  const scalar p,
109  const scalar T
110 ) const
111 {
112  return c0_*pow(T/Tref_, n0_) + EquationOfState::Cv(p, T);
113 }
114 
115 
116 template<class EquationOfState>
118 (
119  const scalar p,
120  const scalar T
121 ) const
122 {
123  return
124  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
125  + EquationOfState::E(p, T);
126 }
127 
128 
129 template<class EquationOfState>
131 (
132  const scalar p,
133  const scalar T
134 ) const
135 {
136  return Es(p, T) + Hf();
137 }
138 
139 
140 template<class EquationOfState>
141 inline Foam::scalar Foam::ePowerThermo<EquationOfState>::Hf() const
142 {
143  return Hf_;
144 }
145 
146 
147 template<class EquationOfState>
149 (
150  const scalar p,
151  const scalar T
152 ) const
153 {
154  return
155  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
156  + EquationOfState::Sv(p, T);
157 }
158 
159 
160 template<class EquationOfState>
162 (
163  const scalar T
164 ) const
165 {
166  return
167  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
168  + Hf() + Pstd/EquationOfState::rho(Pstd, T)
169  - S(Pstd, T)*T;
170 }
171 
172 
173 template<class EquationOfState>
175 (
176  const scalar p,
177  const scalar T
178 ) const
179 {
181  return 0; // + EquationOfState::dCpdT
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
186 
187 template<class EquationOfState>
189 (
191 )
192 {
193  scalar Y1 = this->Y();
194 
195  EquationOfState::operator+=(ct);
196 
197  if (mag(this->Y()) > small)
198  {
199  Y1 /= this->Y();
200  const scalar Y2 = ct.Y()/this->Y();
201 
202  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
203  c0_ = Y1*c0_ + Y2*ct.c0_;
204  n0_ = Y1*n0_ + Y2*ct.n0_;
205  Tref_ = Y1*Tref_ + Y2*ct.Tref_;
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
211 
212 template<class EquationOfState>
213 inline Foam::ePowerThermo<EquationOfState> Foam::operator+
214 (
217 )
218 {
219  EquationOfState eofs
220  (
221  static_cast<const EquationOfState&>(ct1)
222  + static_cast<const EquationOfState&>(ct2)
223  );
224 
225  if (mag(eofs.Y()) < small)
226  {
228  (
229  eofs,
230  ct1.c0_,
231  ct1.n0_,
232  ct1.Tref_,
233  ct1.Hf_
234  );
235  }
236  else
237  {
238  return ePowerThermo<EquationOfState>
239  (
240  eofs,
241  ct1.Y()/eofs.Y()*ct1.c0_
242  + ct2.Y()/eofs.Y()*ct2.c0_,
243  ct1.Y()/eofs.Y()*ct1.n0_
244  + ct2.Y()/eofs.Y()*ct2.n0_,
245  ct1.Y()/eofs.Y()*ct1.Tref_
246  + ct2.Y()/eofs.Y()*ct2.Tref_,
247  ct1.Y()/eofs.Y()*ct1.Hf_
248  + ct2.Y()/eofs.Y()*ct2.Hf_
249  );
250  }
251 }
252 
253 
254 template<class EquationOfState>
255 inline Foam::ePowerThermo<EquationOfState> Foam::operator*
256 (
257  const scalar s,
258  const ePowerThermo<EquationOfState>& ct
259 )
260 {
261  return ePowerThermo<EquationOfState>
262  (
263  s*static_cast<const EquationOfState&>(ct),
264  ct.c0_,
265  ct.n0_,
266  ct.Tref_,
267  ct.Hf_
268  );
269 }
270 
271 
272 template<class EquationOfState>
273 inline Foam::ePowerThermo<EquationOfState> Foam::operator==
274 (
275  const ePowerThermo<EquationOfState>& ct1,
276  const ePowerThermo<EquationOfState>& ct2
277 )
278 {
279  EquationOfState eofs
280  (
281  static_cast<const EquationOfState&>(ct1)
282  == static_cast<const EquationOfState&>(ct2)
283  );
284 
285  return ePowerThermo<EquationOfState>
286  (
287  eofs,
288  ct2.Y()/eofs.Y()*ct2.c0_
289  - ct1.Y()/eofs.Y()*ct1.c0_,
290  ct2.Y()/eofs.Y()*ct2.n0_
291  - ct1.Y()/eofs.Y()*ct1.n0_,
292  ct2.Y()/eofs.Y()*ct2.Tref_
293  - ct1.Y()/eofs.Y()*ct1.Tref_,
294  ct2.Y()/eofs.Y()*ct2.Hf_
295  - ct1.Y()/eofs.Y()*ct1.Hf_
296  );
297 }
298 
299 
300 // ************************************************************************* //
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Internal energy based thermodynamics package using a power function of temperature for the constant h...
Definition: ePowerThermo.H:129
scalar Hf() const
Enthalpy of formation [J/kg].
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: ePowerThermoI.H:97
ePowerThermo(const EquationOfState &st, const scalar c0, const scalar n0, const scalar Tref, const scalar Hf)
Construct from components.
Definition: ePowerThermoI.H:51
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
autoPtr< ePowerThermo > clone() const
Construct and return a clone.
Definition: ePowerThermoI.H:84
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 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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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.
const dimensionedScalar Tstd
Standard temperature.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
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