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-2018 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, const scalar T
120 ) const
121 {
122  return c0_*pow(T/Tref_, n0_) + EquationOfState::Cp(p, T);
123 }
124 
125 
126 template<class EquationOfState>
127 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Ha
128 (
129  const scalar p, const scalar T
130 ) const
131 {
132  return Hs(p, T) + Hc();
133 }
134 
135 
136 template<class EquationOfState>
137 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hs
138 (
139  const scalar p, const scalar T
140 ) const
141 {
142  return
143  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
144  + EquationOfState::H(p, T);
145 }
146 
147 
148 template<class EquationOfState>
149 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hc() const
150 {
151  return Hf_;
152 }
153 
154 
155 template<class EquationOfState>
156 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::S
157 (
158  const scalar p, const scalar T
159 ) const
160 {
161  return
162  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
163  + EquationOfState::S(p, T);
164 }
165 
166 
167 template<class EquationOfState>
169 (
170  const scalar p, const scalar T
171 ) const
172 {
173  // To be implemented
175  return 0;
176 }
177 
178 
179 template<class EquationOfState>
181 (
182  const scalar p, const scalar T
183 ) const
184 {
185  // To be implemented
187  return 0;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
192 
193 template<class EquationOfState>
194 inline void Foam::hPowerThermo<EquationOfState>::operator+=
195 (
197 )
198 {
199  scalar Y1 = this->Y();
200 
201  EquationOfState::operator+=(ct);
202 
203  if (mag(this->Y()) > small)
204  {
205  Y1 /= this->Y();
206  const scalar Y2 = ct.Y()/this->Y();
207 
208  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
209  c0_ = Y1*c0_ + Y2*ct.c0_;
210  n0_ = Y1*n0_ + Y2*ct.n0_;
211  Tref_ = Y1*Tref_ + Y2*ct.Tref_;
212  }
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
217 
218 template<class EquationOfState>
219 inline Foam::hPowerThermo<EquationOfState> Foam::operator+
220 (
223 )
224 {
225  EquationOfState eofs
226  (
227  static_cast<const EquationOfState&>(ct1)
228  + static_cast<const EquationOfState&>(ct2)
229  );
230 
231  if (mag(eofs.Y()) < small)
232  {
234  (
235  eofs,
236  ct1.c0_,
237  ct1.n0_,
238  ct1.Tref_,
239  ct1.Hf_
240  );
241  }
242  else
243  {
245  (
246  eofs,
247  ct1.Y()/eofs.Y()*ct1.c0_
248  + ct2.Y()/eofs.Y()*ct2.c0_,
249  ct1.Y()/eofs.Y()*ct1.n0_
250  + ct2.Y()/eofs.Y()*ct2.n0_,
251  ct1.Y()/eofs.Y()*ct1.Tref_
252  + ct2.Y()/eofs.Y()*ct2.Tref_,
253  ct1.Y()/eofs.Y()*ct1.Hf_
254  + ct2.Y()/eofs.Y()*ct2.Hf_
255  );
256  }
257 }
258 
259 
260 template<class EquationOfState>
261 inline Foam::hPowerThermo<EquationOfState> Foam::operator*
262 (
263  const scalar s,
265 )
266 {
268  (
269  s*static_cast<const EquationOfState&>(ct),
270  ct.c0_,
271  ct.n0_,
272  ct.Tref_,
273  ct.Hf_
274  );
275 }
276 
277 
278 template<class EquationOfState>
279 inline Foam::hPowerThermo<EquationOfState> Foam::operator==
280 (
283 )
284 {
285  EquationOfState eofs
286  (
287  static_cast<const EquationOfState&>(ct1)
288  == static_cast<const EquationOfState&>(ct2)
289  );
290 
292  (
293  eofs,
294  ct2.Y()/eofs.Y()*ct2.c0_
295  - ct1.Y()/eofs.Y()*ct1.c0_,
296  ct2.Y()/eofs.Y()*ct2.n0_
297  - ct1.Y()/eofs.Y()*ct1.n0_,
298  ct2.Y()/eofs.Y()*ct2.Tref_
299  - ct1.Y()/eofs.Y()*ct1.Tref_,
300  ct2.Y()/eofs.Y()*ct2.Hf_
301  - ct1.Y()/eofs.Y()*ct1.Hf_
302  );
303 }
304 
305 
306 // ************************************************************************* //
dictionary dict
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
const dimensionedScalar Tstd
Standard temperature.
const volScalarField & T
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar Hc() const
Chemical enthalpy [J/kg].
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
Heat capacity at constant pressure [J/(kg K)].
Power-function based thermodynamics package templated on EquationOfState.
Definition: hPowerThermo.H:54
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
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:366