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-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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
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 template<class EquationOfState>
48 (
49  const word& name,
50  const hConstThermo& ct
51 )
52 :
53  EquationOfState(name, ct),
54  Cp_(ct.Cp_),
55  Hf_(ct.Hf_),
56  Tref_(ct.Tref_),
57  Hsref_(ct.Hsref_)
58 {}
59 
60 
61 template<class EquationOfState>
64 {
66  (
68  );
69 }
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class EquationOfState>
76 (
77  const scalar T
78 ) const
79 {
80  return T;
81 }
82 
83 
84 template<class EquationOfState>
86 (
87  const scalar p,
88  const scalar T
89 ) const
90 {
91  return Cp_ + EquationOfState::Cp(p, T);
92 }
93 
94 
95 template<class EquationOfState>
97 (
98  const scalar p,
99  const scalar T
100 ) const
101 {
102  return Cp_*(T - Tref_) + Hsref_ + EquationOfState::H(p, T);
103 }
104 
105 
106 template<class EquationOfState>
108 (
109  const scalar p,
110  const scalar T
111 ) const
112 {
113  return Hs(p, T) + Hf();
114 }
115 
116 
117 template<class EquationOfState>
118 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hf() const
119 {
120  return Hf_;
121 }
122 
123 
124 template<class EquationOfState>
126 (
127  const scalar p,
128  const scalar T
129 ) const
130 {
131  return Cp_*log(T/Tstd) + EquationOfState::Sp(p, T);
132 }
133 
134 
135 template<class EquationOfState>
137 (
138  const scalar T
139 ) const
140 {
141  return Cp_*(T - Tref_) + Hsref_ + Hf() - Cp_*T*log(T/Tstd);
142 }
143 
144 
145 template<class EquationOfState>
147 (
148  const scalar p,
149  const scalar T
150 ) const
151 {
152  return 0;
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157 
158 template<class EquationOfState>
160 (
162 )
163 {
164  scalar Y1 = this->Y();
165 
166  EquationOfState::operator+=(ct);
167 
168  if (mag(this->Y()) > small)
169  {
170  if
171  (
173  && notEqual(Tref_, ct.Tref_)
174  )
175  {
177  << "Tref " << Tref_ << " for "
178  << (this->name().size() ? this->name() : "others")
179  << " != " << ct.Tref_ << " for "
180  << (ct.name().size() ? ct.name() : "others")
181  << exit(FatalError);
182  }
183 
184  Y1 /= this->Y();
185  const scalar Y2 = ct.Y()/this->Y();
186 
187  Cp_ = Y1*Cp_ + Y2*ct.Cp_;
188  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
189  Hsref_ = Y1*Hsref_ + Y2*ct.Hsref_;
190  }
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
195 
196 template<class EquationOfState>
197 inline Foam::hConstThermo<EquationOfState> Foam::operator+
198 (
201 )
202 {
203  EquationOfState eofs
204  (
205  static_cast<const EquationOfState&>(ct1)
206  + static_cast<const EquationOfState&>(ct2)
207  );
208 
209  if (mag(eofs.Y()) < small)
210  {
212  (
213  eofs,
214  ct1.Cp_,
215  ct1.Hf_,
216  ct1.Tref_,
217  ct1.Hsref_
218  );
219  }
220  else
221  {
222  if
223  (
224  hConstThermo<EquationOfState>::debug
225  && notEqual(ct1.Tref_, ct2.Tref_)
226  )
227  {
229  << "Tref " << ct1.Tref_ << " for "
230  << (ct1.name().size() ? ct1.name() : "others")
231  << " != " << ct2.Tref_ << " for "
232  << (ct2.name().size() ? ct2.name() : "others")
233  << exit(FatalError);
234  }
235 
236  return hConstThermo<EquationOfState>
237  (
238  eofs,
239  ct1.Y()/eofs.Y()*ct1.Cp_
240  + ct2.Y()/eofs.Y()*ct2.Cp_,
241  ct1.Y()/eofs.Y()*ct1.Hf_
242  + ct2.Y()/eofs.Y()*ct2.Hf_,
243  ct1.Tref_,
244  ct1.Y()/eofs.Y()*ct1.Hsref_
245  + ct2.Y()/eofs.Y()*ct2.Hsref_
246  );
247  }
248 }
249 
250 
251 template<class EquationOfState>
252 inline Foam::hConstThermo<EquationOfState> Foam::operator*
253 (
254  const scalar s,
255  const hConstThermo<EquationOfState>& ct
256 )
257 {
258  return hConstThermo<EquationOfState>
259  (
260  s*static_cast<const EquationOfState&>(ct),
261  ct.Cp_,
262  ct.Hf_,
263  ct.Tref_,
264  ct.Hsref_
265  );
266 }
267 
268 
269 template<class EquationOfState>
270 inline Foam::hConstThermo<EquationOfState> Foam::operator==
271 (
272  const hConstThermo<EquationOfState>& ct1,
273  const hConstThermo<EquationOfState>& ct2
274 )
275 {
276  EquationOfState eofs
277  (
278  static_cast<const EquationOfState&>(ct1)
279  == static_cast<const EquationOfState&>(ct2)
280  );
281 
282  if
283  (
284  hConstThermo<EquationOfState>::debug
285  && notEqual(ct1.Tref_, ct2.Tref_)
286  )
287  {
289  << "Tref " << ct1.Tref_ << " for "
290  << (ct1.name().size() ? ct1.name() : "others")
291  << " != " << ct2.Tref_ << " for "
292  << (ct2.name().size() ? ct2.name() : "others")
293  << exit(FatalError);
294  }
295 
296  return hConstThermo<EquationOfState>
297  (
298  eofs,
299  ct2.Y()/eofs.Y()*ct2.Cp_
300  - ct1.Y()/eofs.Y()*ct1.Cp_,
301  ct2.Y()/eofs.Y()*ct2.Hf_
302  - ct1.Y()/eofs.Y()*ct1.Hf_,
303  ct1.Tref_,
304  ct2.Y()/eofs.Y()*ct2.Hsref_
305  - ct1.Y()/eofs.Y()*ct1.Hsref_
306  );
307 }
308 
309 
310 // ************************************************************************* //
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Enthalpy based thermodynamics package using a constant heat capacity at constant pressure:
Definition: hConstThermo.H:122
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: hConstThermoI.H:76
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hConstThermoI.H:97
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
Definition: hConstThermoI.H:86
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
scalar Ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
autoPtr< hConstThermo > clone() const
Construct and return a clone.
Definition: hConstThermoI.H:63
hConstThermo(const EquationOfState &st, const scalar Cp, const scalar Hf, const scalar Tref, const scalar Hsref)
Construct from components.
Definition: hConstThermoI.H:30
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
A class for handling words, derived from string.
Definition: word.H:62
#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 Tstd
Standard temperature.
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:215
dimensionedScalar log(const dimensionedScalar &ds)
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