hPowerThermoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2016 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 template<class EquationOfState>
107 {
109  (
111  );
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class EquationOfState>
119 (
120  const scalar T
121 ) const
122 {
123  return T;
124 }
125 
126 
127 template<class EquationOfState>
128 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::cp
129 (
130  const scalar p, const scalar T
131 ) const
132 {
133  return c0_*pow(T/Tref_, n0_) + EquationOfState::cp(p, T);
134 }
135 
136 
137 template<class EquationOfState>
138 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::ha
139 (
140  const scalar p, const scalar T
141 ) const
142 {
143  return hs(p, T) + hc();
144 }
145 
146 
147 template<class EquationOfState>
148 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::hs
149 (
150  const scalar p, const scalar T
151 ) const
152 {
153  return
154  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
155  + EquationOfState::h(p, T);
156 }
157 
158 
159 template<class EquationOfState>
160 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::hc() const
161 {
162  return Hf_;
163 }
164 
165 
166 template<class EquationOfState>
167 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::s
168 (
169  const scalar p, const scalar T
170 ) const
171 {
172  return
173  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
174  + EquationOfState::s(p, T);
175 }
176 
177 
178 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
179 
180 template<class EquationOfState>
181 inline void Foam::hPowerThermo<EquationOfState>::operator+=
182 (
184 )
185 {
186  scalar molr1 = this->nMoles();
187 
188  EquationOfState::operator+=(ct);
189  molr1 /= this->nMoles();
190  scalar molr2 = ct.nMoles()/this->nMoles();
191 
192  Hf_ = molr1*Hf_ + molr2*ct.Hf_;
193  c0_ = molr1*c0_ + molr2*ct.c0_;
194  n0_ = molr1*n0_ + molr2*ct.n0_;
195  Tref_ = molr1*Tref_ + molr2*ct.Tref_;
196 }
197 
198 
199 template<class EquationOfState>
200 inline void Foam::hPowerThermo<EquationOfState>::operator-=
201 (
203 )
204 {
205  scalar molr1 = this->nMoles();
206 
207  EquationOfState::operator-=(ct);
208 
209  molr1 /= this->nMoles();
210  scalar molr2 = ct.nMoles()/this->nMoles();
211 
212  Hf_ = molr1*Hf_ - molr2*ct.Hf_;
213  c0_ = (molr1*c0_ - molr2*ct.c0_);
214  n0_ = (molr1*n0_ - molr2*ct.n0_);
215  Tref_ = (molr1*Tref_ - molr2*ct.Tref_);
216 }
217 
218 
219 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
220 
221 template<class EquationOfState>
222 inline Foam::hPowerThermo<EquationOfState> Foam::operator+
223 (
226 )
227 {
228  EquationOfState eofs
229  (
230  static_cast<const EquationOfState&>(ct1)
231  + static_cast<const EquationOfState&>(ct2)
232  );
233 
235  (
236  eofs,
237  ct1.nMoles()/eofs.nMoles()*ct1.c0_
238  + ct2.nMoles()/eofs.nMoles()*ct2.c0_,
239  ct1.nMoles()/eofs.nMoles()*ct1.n0_
240  + ct2.nMoles()/eofs.nMoles()*ct2.n0_,
241  ct1.nMoles()/eofs.nMoles()*ct1.Tref_
242  + ct2.nMoles()/eofs.nMoles()*ct2.Tref_,
243  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
244  + ct2.nMoles()/eofs.nMoles()*ct2.Hf_
245  );
246 }
247 
248 
249 template<class EquationOfState>
250 inline Foam::hPowerThermo<EquationOfState> Foam::operator-
251 (
254 )
255 {
256  EquationOfState eofs
257  (
258  static_cast<const EquationOfState&>(ct1)
259  + static_cast<const EquationOfState&>(ct2)
260  );
261 
263  (
264  eofs,
265  ct1.nMoles()/eofs.nMoles()*ct1.c0_
266  - ct2.nMoles()/eofs.nMoles()*ct2.c0_,
267  ct1.nMoles()/eofs.nMoles()*ct1.n0_
268  - ct2.nMoles()/eofs.nMoles()*ct2.n0_,
269  ct1.nMoles()/eofs.nMoles()*ct1.Tref_
270  - ct2.nMoles()/eofs.nMoles()*ct2.Tref_,
271  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
272  - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
273  );
274 }
275 
276 
277 template<class EquationOfState>
278 inline Foam::hPowerThermo<EquationOfState> Foam::operator*
279 (
280  const scalar s,
282 )
283 {
285  (
286  s*static_cast<const EquationOfState&>(ct),
287  ct.c0_,
288  ct.n0_,
289  ct.Tref_,
290  ct.Hf_
291  );
292 }
293 
294 
295 template<class EquationOfState>
296 inline Foam::hPowerThermo<EquationOfState> Foam::operator==
297 (
300 )
301 {
302  return ct2 - ct1;
303 }
304 
305 
306 // ************************************************************************* //
dictionary dict
autoPtr< hPowerThermo > clone() const
Construct and return a clone.
Definition: hPowerThermoI.H:84
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_.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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 s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
A class for handling words, derived from string.
Definition: word.H:59
static autoPtr< hPowerThermo > New(Istream &is)
Selector from Istream.
Definition: hPowerThermoI.H:95
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionedScalar Tstd
Standard temperature.
const volScalarField & cp
const volScalarField & T
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kmol].
Power-function based thermodynamics package templated on EquationOfState.
Definition: hPowerThermo.H:54
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
scalar hc() const
Chemical enthalpy [J/kg].