eConstThermoI.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-2020 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 template<class EquationOfState>
30 (
31  const EquationOfState& st,
32  const scalar Cv,
33  const scalar Hf,
34  const scalar Tref,
35  const scalar Esref
36 )
37 :
38  EquationOfState(st),
39  Cv_(Cv),
40  Hf_(Hf),
41  Tref_(Tref),
42  Esref_(Esref)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class EquationOfState>
50 (
51  const word& name,
52  const eConstThermo& ct
53 )
54 :
55  EquationOfState(name, ct),
56  Cv_(ct.Cv_),
57  Hf_(ct.Hf_),
58  Tref_(ct.Tref_),
59  Esref_(ct.Esref_)
60 {}
61 
62 
63 template<class EquationOfState>
66 {
68  (
70  );
71 }
72 
73 
74 template<class EquationOfState>
77 {
79  (
81  );
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class EquationOfState>
89 (
90  const scalar T
91 ) const
92 {
93  return T;
94 }
95 
96 
97 template<class EquationOfState>
99 (
100  const scalar p,
101  const scalar T
102 ) const
103 {
104  return Cv_ + EquationOfState::Cv(p, T);
105 }
106 
107 
108 template<class EquationOfState>
109 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Es
110 (
111  const scalar p,
112  const scalar T
113 ) const
114 {
115  return Cv_*(T - Tref_) + Esref_ + EquationOfState::E(p, T);
116 }
117 
118 
119 template<class EquationOfState>
120 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Hf() const
121 {
122  return Hf_;
123 }
124 
125 
126 template<class EquationOfState>
127 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Ea
128 (
129  const scalar p,
130  const scalar T
131 ) const
132 {
133  return Es(p, T) + Hf();
134 }
135 
136 
137 template<class EquationOfState>
138 inline Foam::scalar Foam::eConstThermo<EquationOfState>::S
139 (
140  const scalar p,
141  const scalar T
142 ) const
143 {
144  return Cp(p, T)*log(T/Tstd) + EquationOfState::Sv(p, T);
145 }
146 
147 
148 template<class EquationOfState>
150 (
151  const scalar T
152 ) const
153 {
154  return
155  Cv_*(T - Tref_) + Esref_ + Hf() + Pstd/EquationOfState::rho(Pstd, T)
156  - S(Pstd, T)*T;
157 }
158 
159 
160 template<class EquationOfState>
162 (
163  const scalar p,
164  const scalar T
165 ) const
166 {
168  return 0; // EquationOfState::dCpdT
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
173 
174 template<class EquationOfState>
175 inline void Foam::eConstThermo<EquationOfState>::operator+=
176 (
178 )
179 {
180  scalar Y1 = this->Y();
181 
182  EquationOfState::operator+=(ct);
183 
184  if (mag(this->Y()) > small)
185  {
186  if
187  (
189  && notEqual(Tref_, ct.Tref_)
190  )
191  {
193  << "Tref " << Tref_ << " for "
194  << (this->name().size() ? this->name() : "others")
195  << " != " << ct.Tref_ << " for "
196  << (ct.name().size() ? ct.name() : "others")
197  << exit(FatalError);
198  }
199 
200  Y1 /= this->Y();
201  const scalar Y2 = ct.Y()/this->Y();
202 
203  Cv_ = Y1*Cv_ + Y2*ct.Cv_;
204  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
205  Esref_ = Y1*Esref_ + Y2*ct.Esref_;
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
211 
212 template<class EquationOfState>
213 inline Foam::eConstThermo<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.Cv_,
231  ct1.Hf_,
232  ct1.Tref_,
233  ct1.Esref_
234  );
235  }
236  else
237  {
238  if
239  (
241  && notEqual(ct1.Tref_, ct2.Tref_)
242  )
243  {
245  << "Tref " << ct1.Tref_ << " for "
246  << (ct1.name().size() ? ct1.name() : "others")
247  << " != " << ct2.Tref_ << " for "
248  << (ct2.name().size() ? ct2.name() : "others")
249  << exit(FatalError);
250  }
251 
253  (
254  eofs,
255  ct1.Y()/eofs.Y()*ct1.Cv_
256  + ct2.Y()/eofs.Y()*ct2.Cv_,
257  ct1.Y()/eofs.Y()*ct1.Hf_
258  + ct2.Y()/eofs.Y()*ct2.Hf_,
259  ct1.Tref_,
260  ct1.Y()/eofs.Y()*ct1.Esref_
261  + ct2.Y()/eofs.Y()*ct2.Esref_
262  );
263  }
264 }
265 
266 
267 template<class EquationOfState>
268 inline Foam::eConstThermo<EquationOfState> Foam::operator*
269 (
270  const scalar s,
272 )
273 {
275  (
276  s*static_cast<const EquationOfState&>(ct),
277  ct.Cv_,
278  ct.Hf_,
279  ct.Tref_,
280  ct.Esref_
281  );
282 }
283 
284 
285 template<class EquationOfState>
286 inline Foam::eConstThermo<EquationOfState> Foam::operator==
287 (
290 )
291 {
292  EquationOfState eofs
293  (
294  static_cast<const EquationOfState&>(ct1)
295  == static_cast<const EquationOfState&>(ct2)
296  );
297 
298  if
299  (
301  && notEqual(ct1.Tref_, ct2.Tref_)
302  )
303  {
305  << "Tref " << ct1.Tref_ << " for "
306  << (ct1.name().size() ? ct1.name() : "others")
307  << " != " << ct2.Tref_ << " for "
308  << (ct2.name().size() ? ct2.name() : "others")
309  << exit(FatalError);
310  }
311 
313  (
314  eofs,
315  ct2.Y()/eofs.Y()*ct2.Cv_
316  - ct1.Y()/eofs.Y()*ct1.Cv_,
317  ct2.Y()/eofs.Y()*ct2.Hf_
318  - ct1.Y()/eofs.Y()*ct1.Hf_,
319  ct1.Tref_,
320  ct2.Y()/eofs.Y()*ct2.Esref_
321  - ct1.Y()/eofs.Y()*ct1.Esref_
322  );
323 }
324 
325 
326 // ************************************************************************* //
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: eConstThermoI.H:89
dictionary dict
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:215
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const dimensionedScalar & Tstd
Standard temperature.
const volScalarField & Cv
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
scalar Hf() const
Enthalpy of formation [J/kg].
Constant properties thermodynamics package templated on an equation of state.
Definition: eConstThermo.H:46
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [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
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< eConstThermo > clone() const
Construct and return a clone.
Definition: eConstThermoI.H:65
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/kg/K].
Definition: eConstThermoI.H:99
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar & Pstd
Standard pressure.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static autoPtr< eConstThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: eConstThermoI.H:76
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366