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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "hPowerThermo.H"
29 #include "specie.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class EquationOfState>
35 (
36  const scalar T
37 ) const
38 {
39  if (T < 0)
40  {
42  << "attempt to evaluate hPowerThermo<EquationOfState>"
43  " for negative temperature " << T
44  << abort(FatalError);
45  }
46 }
47 
48 
49 template<class EquationOfState>
51 (
52  const word& name,
53  const hPowerThermo& jt
54 )
55 :
56  EquationOfState(name, jt),
57  c0_(jt.c0_),
58  n0_(jt.n0_),
59  Tref_(jt.Tref_),
60  Hf_(jt.Hf_)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
66 template<class EquationOfState>
68 (
69  const EquationOfState& st,
70  const scalar c0,
71  const scalar n0,
72  const scalar Tref,
73  const scalar Hf
74 )
75 :
76  EquationOfState(st),
77  c0_(c0),
78  n0_(n0),
79  Tref_(Tref),
80  Hf_(Hf)
81 {}
82 
83 
84 template<class EquationOfState>
87 {
89 }
90 
91 
92 template<class EquationOfState>
95 {
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 template<class EquationOfState>
104 (
105  const scalar T
106 ) const
107 {
108  return T;
109 }
110 
111 
112 template<class EquationOfState>
113 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Cp
114 (
115  const scalar p, const scalar T
116 ) const
117 {
118  return c0_*pow(T/Tref_, n0_) + EquationOfState::Cp(p, T);
119 }
120 
121 
122 template<class EquationOfState>
123 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Ha
124 (
125  const scalar p, const scalar T
126 ) const
127 {
128  return Hs(p, T) + Hc();
129 }
130 
131 
132 template<class EquationOfState>
133 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hs
134 (
135  const scalar p, const scalar T
136 ) const
137 {
138  return
139  c0_*(pow(T, n0_ + 1) - pow(Tstd, n0_ + 1))/(pow(Tref_, n0_)*(n0_ + 1))
140  + EquationOfState::H(p, T);
141 }
142 
143 
144 template<class EquationOfState>
145 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::Hc() const
146 {
147  return Hf_;
148 }
149 
150 
151 template<class EquationOfState>
152 inline Foam::scalar Foam::hPowerThermo<EquationOfState>::S
153 (
154  const scalar p, const scalar T
155 ) const
156 {
157  return
158  c0_*(pow(T, n0_) - pow(Tstd, n0_))/(pow(Tref_, n0_)*n0_)
159  + EquationOfState::S(p, T);
160 }
161 
162 
163 template<class EquationOfState>
165 (
166  const scalar p, const scalar T
167 ) const
168 {
169  // To be implemented
171  return 0;
172 }
173 
174 
175 template<class EquationOfState>
177 (
178  const scalar p, const scalar T
179 ) const
180 {
181  // To be implemented
183  return 0;
184 }
185 
186 
187 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
188 
189 template<class EquationOfState>
191 (
193 )
194 {
195  scalar Y1 = this->Y();
196 
197  EquationOfState::operator+=(ct);
198 
199  if (mag(this->Y()) > SMALL)
200  {
201  Y1 /= this->Y();
202  const scalar Y2 = ct.Y()/this->Y();
203 
204  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
205  c0_ = Y1*c0_ + Y2*ct.c0_;
206  n0_ = Y1*n0_ + Y2*ct.n0_;
207  Tref_ = Y1*Tref_ + Y2*ct.Tref_;
208  }
209 }
210 
211 
212 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
213 
214 template<class EquationOfState>
215 inline Foam::hPowerThermo<EquationOfState> Foam::operator+
216 (
219 )
220 {
221  EquationOfState eofs
222  (
223  static_cast<const EquationOfState&>(ct1)
224  + static_cast<const EquationOfState&>(ct2)
225  );
226 
227  if (mag(eofs.Y()) < SMALL)
228  {
230  (
231  eofs,
232  ct1.c0_,
233  ct1.n0_,
234  ct1.Tref_,
235  ct1.Hf_
236  );
237  }
238  else
239  {
240  return hPowerThermo<EquationOfState>
241  (
242  eofs,
243  ct1.Y()/eofs.Y()*ct1.c0_
244  + ct2.Y()/eofs.Y()*ct2.c0_,
245  ct1.Y()/eofs.Y()*ct1.n0_
246  + ct2.Y()/eofs.Y()*ct2.n0_,
247  ct1.Y()/eofs.Y()*ct1.Tref_
248  + ct2.Y()/eofs.Y()*ct2.Tref_,
249  ct1.Y()/eofs.Y()*ct1.Hf_
250  + ct2.Y()/eofs.Y()*ct2.Hf_
251  );
252  }
253 }
254 
255 
256 template<class EquationOfState>
257 inline Foam::hPowerThermo<EquationOfState> Foam::operator*
258 (
259  const scalar s,
260  const hPowerThermo<EquationOfState>& ct
261 )
262 {
263  return hPowerThermo<EquationOfState>
264  (
265  s*static_cast<const EquationOfState&>(ct),
266  ct.c0_,
267  ct.n0_,
268  ct.Tref_,
269  ct.Hf_
270  );
271 }
272 
273 
274 template<class EquationOfState>
275 inline Foam::hPowerThermo<EquationOfState> Foam::operator==
276 (
277  const hPowerThermo<EquationOfState>& ct1,
278  const hPowerThermo<EquationOfState>& ct2
279 )
280 {
281  EquationOfState eofs
282  (
283  static_cast<const EquationOfState&>(ct1)
284  == static_cast<const EquationOfState&>(ct2)
285  );
286 
287  return hPowerThermo<EquationOfState>
288  (
289  eofs,
290  ct2.Y()/eofs.Y()*ct2.c0_
291  - ct1.Y()/eofs.Y()*ct1.c0_,
292  ct2.Y()/eofs.Y()*ct2.n0_
293  - ct1.Y()/eofs.Y()*ct1.n0_,
294  ct2.Y()/eofs.Y()*ct2.Tref_
295  - ct1.Y()/eofs.Y()*ct1.Tref_,
296  ct2.Y()/eofs.Y()*ct2.Hf_
297  - ct1.Y()/eofs.Y()*ct1.Hf_
298  );
299 }
300 
301 
302 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::hPowerThermo
Power-function based thermodynamics package templated on EquationOfState.
Definition: hPowerThermo.H:59
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
hPowerThermo.H
Hs
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:17
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::hPowerThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hPowerThermoI.H:104
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
specie.H
Foam::hPowerThermo::Hs
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hPowerThermoI.H:134
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:419
Foam::hPowerThermo::clone
autoPtr< hPowerThermo > clone() const
Construct and return a clone.
Definition: hPowerThermoI.H:86
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::hPowerThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: hPowerThermoI.H:145
Foam::hPowerThermo::Ha
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
Definition: hPowerThermoI.H:124
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::hPowerThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: hPowerThermoI.H:153
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::hPowerThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: hPowerThermoI.H:177
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::hPowerThermo::New
static autoPtr< hPowerThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: hPowerThermoI.H:94
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::hPowerThermo::Cp
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: hPowerThermoI.H:114
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::hPowerThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: hPowerThermoI.H:165