All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
hConstThermoI.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 Cp,
33  const scalar Hf,
34  const scalar Tref,
35  const scalar Href
36 )
37 :
38  EquationOfState(st),
39  Cp_(Cp),
40  Hf_(Hf),
41  Tref_(Tref),
42  Hsref_(Href)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class EquationOfState>
50 (
51  const word& name,
52  const hConstThermo& ct
53 )
54 :
55  EquationOfState(name, ct),
56  Cp_(ct.Cp_),
57  Hf_(ct.Hf_),
58  Tref_(ct.Tref_),
59  Hsref_(ct.Hsref_)
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 Cp_ + EquationOfState::Cp(p, T);
105 }
106 
107 
108 template<class EquationOfState>
109 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hs
110 (
111  const scalar p,
112  const scalar T
113 ) const
114 {
115  return Cp_*(T - Tref_) + Hsref_ + EquationOfState::H(p, T);
116 }
117 
118 
119 template<class EquationOfState>
120 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Ha
121 (
122  const scalar p,
123  const scalar T
124 ) const
125 {
126  return Hs(p, T) + Hf();
127 }
128 
129 
130 template<class EquationOfState>
131 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hf() const
132 {
133  return Hf_;
134 }
135 
136 
137 template<class EquationOfState>
138 inline Foam::scalar Foam::hConstThermo<EquationOfState>::S
139 (
140  const scalar p,
141  const scalar T
142 ) const
143 {
144  return Cp_*log(T/Tstd) + EquationOfState::Sp(p, T);
145 }
146 
147 
148 template<class EquationOfState>
150 (
151  const scalar T
152 ) const
153 {
154  return Cp_*(T - Tref_) + Hsref_ + Hf() - Cp_*T*log(T/Tstd);
155 }
156 
157 
158 template<class EquationOfState>
160 (
161  const scalar p,
162  const scalar T
163 ) const
164 {
165  return 0;
166 }
167 
168 
169 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
170 
171 template<class EquationOfState>
172 inline void Foam::hConstThermo<EquationOfState>::operator+=
173 (
175 )
176 {
177  scalar Y1 = this->Y();
178 
179  EquationOfState::operator+=(ct);
180 
181  if (mag(this->Y()) > small)
182  {
183  if
184  (
186  && notEqual(Tref_, ct.Tref_)
187  )
188  {
190  << "Tref " << Tref_ << " for "
191  << (this->name().size() ? this->name() : "others")
192  << " != " << ct.Tref_ << " for "
193  << (ct.name().size() ? ct.name() : "others")
194  << exit(FatalError);
195  }
196 
197  Y1 /= this->Y();
198  const scalar Y2 = ct.Y()/this->Y();
199 
200  Cp_ = Y1*Cp_ + Y2*ct.Cp_;
201  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
202  Hsref_ = Y1*Hsref_ + Y2*ct.Hsref_;
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
208 
209 template<class EquationOfState>
210 inline Foam::hConstThermo<EquationOfState> Foam::operator+
211 (
214 )
215 {
216  EquationOfState eofs
217  (
218  static_cast<const EquationOfState&>(ct1)
219  + static_cast<const EquationOfState&>(ct2)
220  );
221 
222  if (mag(eofs.Y()) < small)
223  {
225  (
226  eofs,
227  ct1.Cp_,
228  ct1.Hf_,
229  ct1.Tref_,
230  ct1.Hsref_
231  );
232  }
233  else
234  {
235  if
236  (
238  && notEqual(ct1.Tref_, ct2.Tref_)
239  )
240  {
242  << "Tref " << ct1.Tref_ << " for "
243  << (ct1.name().size() ? ct1.name() : "others")
244  << " != " << ct2.Tref_ << " for "
245  << (ct2.name().size() ? ct2.name() : "others")
246  << exit(FatalError);
247  }
248 
250  (
251  eofs,
252  ct1.Y()/eofs.Y()*ct1.Cp_
253  + ct2.Y()/eofs.Y()*ct2.Cp_,
254  ct1.Y()/eofs.Y()*ct1.Hf_
255  + ct2.Y()/eofs.Y()*ct2.Hf_,
256  ct1.Tref_,
257  ct1.Y()/eofs.Y()*ct1.Hsref_
258  + ct2.Y()/eofs.Y()*ct2.Hsref_
259  );
260  }
261 }
262 
263 
264 template<class EquationOfState>
265 inline Foam::hConstThermo<EquationOfState> Foam::operator*
266 (
267  const scalar s,
269 )
270 {
272  (
273  s*static_cast<const EquationOfState&>(ct),
274  ct.Cp_,
275  ct.Hf_,
276  ct.Tref_,
277  ct.Hsref_
278  );
279 }
280 
281 
282 template<class EquationOfState>
283 inline Foam::hConstThermo<EquationOfState> Foam::operator==
284 (
287 )
288 {
289  EquationOfState eofs
290  (
291  static_cast<const EquationOfState&>(ct1)
292  == static_cast<const EquationOfState&>(ct2)
293  );
294 
295  if
296  (
298  && notEqual(ct1.Tref_, ct2.Tref_)
299  )
300  {
302  << "Tref " << ct1.Tref_ << " for "
303  << (ct1.name().size() ? ct1.name() : "others")
304  << " != " << ct2.Tref_ << " for "
305  << (ct2.name().size() ? ct2.name() : "others")
306  << exit(FatalError);
307  }
308 
310  (
311  eofs,
312  ct2.Y()/eofs.Y()*ct2.Cp_
313  - ct1.Y()/eofs.Y()*ct1.Cp_,
314  ct2.Y()/eofs.Y()*ct2.Hf_
315  - ct1.Y()/eofs.Y()*ct1.Hf_,
316  ct1.Tref_,
317  ct2.Y()/eofs.Y()*ct2.Hsref_
318  - ct1.Y()/eofs.Y()*ct1.Hsref_
319  );
320 }
321 
322 
323 // ************************************************************************* //
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
dictionary dict
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:215
Enthalpy based thermodynamics package using a constant heat capacity at constant pressure: ...
Definition: hConstThermo.H:83
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hConstThermoI.H:89
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:156
const dimensionedScalar Tstd
Standard temperature.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
autoPtr< hConstThermo > clone() const
Construct and return a clone.
Definition: hConstThermoI.H:65
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: hConstThermoI.H:76
scalar Hf() const
Enthalpy of formation [J/kg].
scalar Ha(const scalar p, const scalar T) const
Absolute enthalpy [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 Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
Definition: hConstThermoI.H:99
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)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52