All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
hPowerThermoI.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) 2012-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 #include "hPowerThermo.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 hPowerThermo<EquationOfState>"
41  " for negative temperature " << T
42  << abort(FatalError);
43  }
44 }
45 
46 
47 template<class EquationOfState>
49 (
50  const word& name,
51  const hPowerThermo& jt
52 )
53 :
54  EquationOfState(name, jt),
55  c0_(jt.c0_),
56  n0_(jt.n0_),
57  Tref_(jt.Tref_),
58  Hf_(jt.Hf_)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
64 template<class EquationOfState>
66 (
67  const EquationOfState& st,
68  const scalar c0,
69  const scalar n0,
70  const scalar Tref,
71  const scalar Hf
72 )
73 :
74  EquationOfState(st),
75  c0_(c0),
76  n0_(n0),
77  Tref_(Tref),
78  Hf_(Hf)
79 {}
80 
81 
82 template<class EquationOfState>
85 {
87  (
89  );
90 }
91 
92 
93 template<class EquationOfState>
96 {
98  (
100  );
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class EquationOfState>
108 (
109  const scalar T
110 ) const
111 {
112  return T;
113 }
114 
115 
116 template<class EquationOfState>
117 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Cp
118 (
119  const scalar p,
120  const scalar T
121 ) const
122 {
123  return c0_*pow(T/Tref_, n0_) + EquationOfState::Cp(p, T);
124 }
125 
126 
127 template<class EquationOfState>
128 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hs
129 (
130  const scalar p,
131  const scalar T
132 ) const
133 {
134  return
135  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
136  + EquationOfState::H(p, T);
137 }
138 
139 
140 template<class EquationOfState>
141 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Ha
142 (
143  const scalar p,
144  const scalar T
145 ) const
146 {
147  return Hs(p, T) + Hf();
148 }
149 
150 
151 template<class EquationOfState>
152 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hf() const
153 {
154  return Hf_;
155 }
156 
157 
158 template<class EquationOfState>
159 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::S
160 (
161  const scalar p,
162  const scalar T
163 ) const
164 {
165  return
166  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
167  + EquationOfState::Sp(p, T);
168 }
169 
170 
171 template<class EquationOfState>
173 (
174  const scalar T
175 ) const
176 {
177  return
178  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
179  + Hf()
180  - c0_*(pow(T, n0_) - pow(Tstd, n0_))*T/(pow(Tref_, n0_)*n0_);
181 }
182 
183 
184 template<class EquationOfState>
186 (
187  const scalar p,
188  const scalar T
189 ) const
190 {
191  // To be implemented
193  return 0;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
198 
199 template<class EquationOfState>
200 inline void Foam::hPowerThermo<EquationOfState>::operator+=
201 (
203 )
204 {
205  scalar Y1 = this->Y();
206 
207  EquationOfState::operator+=(ct);
208 
209  if (mag(this->Y()) > small)
210  {
211  Y1 /= this->Y();
212  const scalar Y2 = ct.Y()/this->Y();
213 
214  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
215  c0_ = Y1*c0_ + Y2*ct.c0_;
216  n0_ = Y1*n0_ + Y2*ct.n0_;
217  Tref_ = Y1*Tref_ + Y2*ct.Tref_;
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
223 
224 template<class EquationOfState>
225 inline Foam::hPowerThermo<EquationOfState> Foam::operator+
226 (
229 )
230 {
231  EquationOfState eofs
232  (
233  static_cast<const EquationOfState&>(ct1)
234  + static_cast<const EquationOfState&>(ct2)
235  );
236 
237  if (mag(eofs.Y()) < small)
238  {
240  (
241  eofs,
242  ct1.c0_,
243  ct1.n0_,
244  ct1.Tref_,
245  ct1.Hf_
246  );
247  }
248  else
249  {
251  (
252  eofs,
253  ct1.Y()/eofs.Y()*ct1.c0_
254  + ct2.Y()/eofs.Y()*ct2.c0_,
255  ct1.Y()/eofs.Y()*ct1.n0_
256  + ct2.Y()/eofs.Y()*ct2.n0_,
257  ct1.Y()/eofs.Y()*ct1.Tref_
258  + ct2.Y()/eofs.Y()*ct2.Tref_,
259  ct1.Y()/eofs.Y()*ct1.Hf_
260  + ct2.Y()/eofs.Y()*ct2.Hf_
261  );
262  }
263 }
264 
265 
266 template<class EquationOfState>
267 inline Foam::hPowerThermo<EquationOfState> Foam::operator*
268 (
269  const scalar s,
271 )
272 {
274  (
275  s*static_cast<const EquationOfState&>(ct),
276  ct.c0_,
277  ct.n0_,
278  ct.Tref_,
279  ct.Hf_
280  );
281 }
282 
283 
284 template<class EquationOfState>
285 inline Foam::hPowerThermo<EquationOfState> Foam::operator==
286 (
289 )
290 {
291  EquationOfState eofs
292  (
293  static_cast<const EquationOfState&>(ct1)
294  == static_cast<const EquationOfState&>(ct2)
295  );
296 
298  (
299  eofs,
300  ct2.Y()/eofs.Y()*ct2.c0_
301  - ct1.Y()/eofs.Y()*ct1.c0_,
302  ct2.Y()/eofs.Y()*ct2.n0_
303  - ct1.Y()/eofs.Y()*ct1.n0_,
304  ct2.Y()/eofs.Y()*ct2.Tref_
305  - ct1.Y()/eofs.Y()*ct1.Tref_,
306  ct2.Y()/eofs.Y()*ct2.Hf_
307  - ct1.Y()/eofs.Y()*ct1.Hf_
308  );
309 }
310 
311 
312 // ************************************************************************* //
dictionary dict
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
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
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar Ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
scalar S(const scalar p, const scalar T) const
Entropy [J/kg/K].
autoPtr< hPowerThermo > clone() const
Construct and return a clone.
Definition: hPowerThermoI.H:84
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))
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
A class for handling words, derived from string.
Definition: word.H:59
scalar Hf() const
Enthalpy of formation [J/kg].
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const volScalarField & T
static autoPtr< hPowerThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: hPowerThermoI.H:95
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353