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-2015 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_);
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 }
156 
157 
158 template<class EquationOfState>
159 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::hc() const
160 {
161  return Hf_;
162 }
163 
164 
165 template<class EquationOfState>
166 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::s
167 (
168  const scalar p, const scalar T
169 ) const
170 {
171  return
172  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
173  + EquationOfState::s(p, T);
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
178 
179 template<class EquationOfState>
181 (
183 )
184 {
185  scalar molr1 = this->nMoles();
186 
187  EquationOfState::operator+=(ct);
188  molr1 /= this->nMoles();
189  scalar molr2 = ct.nMoles()/this->nMoles();
190 
191  Hf_ = molr1*Hf_ + molr2*ct.Hf_;
192  c0_ = molr1*c0_ + molr2*ct.c0_;
193  n0_ = molr1*n0_ + molr2*ct.n0_;
194  Tref_ = molr1*Tref_ + molr2*ct.Tref_;
195 }
196 
197 
198 template<class EquationOfState>
200 (
202 )
203 {
204  scalar molr1 = this->nMoles();
205 
206  EquationOfState::operator-=(ct);
207 
208  molr1 /= this->nMoles();
209  scalar molr2 = ct.nMoles()/this->nMoles();
210 
211  Hf_ = molr1*Hf_ - molr2*ct.Hf_;
212  c0_ = (molr1*c0_ - molr2*ct.c0_);
213  n0_ = (molr1*n0_ - molr2*ct.n0_);
214  Tref_ = (molr1*Tref_ - molr2*ct.Tref_);
215 }
216 
217 
218 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
219 
220 template<class EquationOfState>
221 inline Foam::hPowerThermo<EquationOfState> Foam::operator+
222 (
225 )
226 {
227  EquationOfState eofs
228  (
229  static_cast<const EquationOfState&>(ct1)
230  + static_cast<const EquationOfState&>(ct2)
231  );
232 
234  (
235  eofs,
236  ct1.nMoles()/eofs.nMoles()*ct1.c0_
237  + ct2.nMoles()/eofs.nMoles()*ct2.c0_,
238  ct1.nMoles()/eofs.nMoles()*ct1.n0_
239  + ct2.nMoles()/eofs.nMoles()*ct2.n0_,
240  ct1.nMoles()/eofs.nMoles()*ct1.Tref_
241  + ct2.nMoles()/eofs.nMoles()*ct2.Tref_,
242  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
243  + ct2.nMoles()/eofs.nMoles()*ct2.Hf_
244  );
245 }
246 
247 
248 template<class EquationOfState>
249 inline Foam::hPowerThermo<EquationOfState> Foam::operator-
250 (
251  const hPowerThermo<EquationOfState>& ct1,
252  const hPowerThermo<EquationOfState>& ct2
253 )
254 {
255  EquationOfState eofs
256  (
257  static_cast<const EquationOfState&>(ct1)
258  + static_cast<const EquationOfState&>(ct2)
259  );
260 
261  return hPowerThermo<EquationOfState>
262  (
263  eofs,
264  ct1.nMoles()/eofs.nMoles()*ct1.c0_
265  - ct2.nMoles()/eofs.nMoles()*ct2.c0_,
266  ct1.nMoles()/eofs.nMoles()*ct1.n0_
267  - ct2.nMoles()/eofs.nMoles()*ct2.n0_,
268  ct1.nMoles()/eofs.nMoles()*ct1.Tref_
269  - ct2.nMoles()/eofs.nMoles()*ct2.Tref_,
270  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
271  - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
272  );
273 }
274 
275 
276 template<class EquationOfState>
277 inline Foam::hPowerThermo<EquationOfState> Foam::operator*
278 (
279  const scalar s,
280  const hPowerThermo<EquationOfState>& ct
281 )
282 {
283  return hPowerThermo<EquationOfState>
284  (
285  s*static_cast<const EquationOfState&>(ct),
286  ct.c0_,
287  ct.n0_,
288  ct.Tref_,
289  ct.Hf_
290  );
291 }
292 
293 
294 template<class EquationOfState>
295 inline Foam::hPowerThermo<EquationOfState> Foam::operator==
296 (
297  const hPowerThermo<EquationOfState>& ct1,
298  const hPowerThermo<EquationOfState>& ct2
299 )
300 {
301  return ct2 - ct1;
302 }
303 
304 
305 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::hPowerThermo
Power-function based thermodynamics package templated on EquationOfState.
Definition: hPowerThermo.H:54
hPowerThermo.H
Foam::hPowerThermo::hPowerThermo
hPowerThermo(const EquationOfState &st, const scalar c0, const scalar n0, const scalar Tref, const scalar Hf)
Construct from components.
Definition: hPowerThermoI.H:66
Foam::hPowerThermo::n0_
scalar n0_
Definition: hPowerThermo.H:106
Foam::hPowerThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hPowerThermoI.H:119
Foam::hPowerThermo::New
static autoPtr< hPowerThermo > New(Istream &is)
Selector from Istream.
Definition: hPowerThermoI.H:95
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:47
specie.H
Foam::hPowerThermo::ha
scalar ha(const scalar p, const scalar T) const
Absolute enthalpy [J/kmol].
Definition: hPowerThermoI.H:139
Foam::hPowerThermo::checkT
void checkT(const scalar T) const
Check given temperature is within the range of the fitted coeffs.
Definition: hPowerThermoI.H:33
Foam::hPowerThermo::cp
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: hPowerThermoI.H:129
Foam::hPowerThermo::hs
scalar hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hPowerThermoI.H:149
Foam::hPowerThermo::clone
autoPtr< hPowerThermo > clone() const
Construct and return a clone.
Definition: hPowerThermoI.H:84
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::hPowerThermo::c0_
scalar c0_
Definition: hPowerThermo.H:105
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::hPowerThermo::Hf_
scalar Hf_
Definition: hPowerThermo.H:108
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::hPowerThermo::hc
scalar hc() const
Chemical enthalpy [J/kg].
Definition: hPowerThermoI.H:159
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
s
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))
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
T
const volScalarField & T
Definition: createFields.H:25
Foam::hPowerThermo::s
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
Definition: hPowerThermoI.H:167
Foam::hPowerThermo::Tref_
scalar Tref_
Definition: hPowerThermo.H:107
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47