Loading...
Searching...
No Matches
hConstThermoI.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29
30template<class EquationOfState>
31inline Foam::hConstThermo<EquationOfState>::hConstThermo
32(
33 const EquationOfState& st,
34 const scalar cp,
35 const scalar hf,
36 const scalar Tref,
37 const scalar Href
38)
39:
40 EquationOfState(st),
41 Cp_(cp),
42 Hf_(hf),
43 Tref_(Tref),
44 Hsref_(Href)
45{}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
50template<class EquationOfState>
51inline Foam::hConstThermo<EquationOfState>::hConstThermo
52(
53 const word& name,
54 const hConstThermo& ct
55)
56:
57 EquationOfState(name, ct),
58 Cp_(ct.Cp_),
59 Hf_(ct.Hf_),
60 Tref_(ct.Tref_),
61 Hsref_(ct.Hsref_)
62{}
63
64
65template<class EquationOfState>
72
73template<class EquationOfState>
76{
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
83template<class EquationOfState>
85(
86 const scalar T
87) const
88{
89 return T;
90}
91
92
93template<class EquationOfState>
95(
96 const scalar p,
97 const scalar T
98) const
99{
100 return Cp_ + EquationOfState::Cp(p, T);
101}
102
103
104template<class EquationOfState>
106(
107 const scalar p, const scalar T
108) const
109{
110 return Hs(p, T) + Hc();
111}
112
113
114template<class EquationOfState>
116(
117 const scalar p, const scalar T
118) const
119{
120 return Cp_*(T - Tref_) + Hsref_ + EquationOfState::H(p, T);
121}
122
123
124template<class EquationOfState>
126{
127 return Hf_;
128}
129
130
131template<class EquationOfState>
133(
134 const scalar p, const scalar T
135) const
136{
137 return Cp_*log(T/Tstd) + EquationOfState::S(p, T);
138}
139
140
141template<class EquationOfState>
143(
144 const scalar T
145) const
146{
147 return Cp_*(T - Tref_) + Hsref_ + Hc() - Cp_*T*log(T/Tstd);
148}
149
150
151template<class EquationOfState>
153(
154 const scalar p, const scalar T
155) const
157 return 0;
158}
159
160// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
161
162template<class EquationOfState>
164(
166)
167{
168 scalar Y1 = this->Y();
169
170 EquationOfState::operator+=(ct);
171
172 if (mag(this->Y()) > SMALL)
173 {
174 Y1 /= this->Y();
175 scalar Y2 = ct.Y()/this->Y();
176
177 Cp_ = Y1*Cp_ + Y2*ct.Cp_;
178 Hf_ = Y1*Hf_ + Y2*ct.Hf_;
179 Hsref_ = Y1*Hsref_ + Y2*ct.Hsref_;
180 }
181}
182
183
184// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
185
186template<class EquationOfState>
187inline Foam::hConstThermo<EquationOfState> Foam::operator+
188(
189 const hConstThermo<EquationOfState>& ct1,
190 const hConstThermo<EquationOfState>& ct2
191)
192{
193 EquationOfState eofs
194 (
195 static_cast<const EquationOfState&>(ct1)
196 + static_cast<const EquationOfState&>(ct2)
197 );
198
199 if (mag(eofs.Y()) < SMALL)
200 {
201 return hConstThermo<EquationOfState>
202 (
203 eofs,
204 ct1.Cp_,
205 ct1.Hf_,
206 ct1.Tref_,
207 ct1.Hsref_
208 );
209 }
210 else
211 {
212 return hConstThermo<EquationOfState>
213 (
214 eofs,
215 ct1.Y()/eofs.Y()*ct1.Cp_
216 + ct2.Y()/eofs.Y()*ct2.Cp_,
217 ct1.Y()/eofs.Y()*ct1.Hf_
218 + ct2.Y()/eofs.Y()*ct2.Hf_,
219 ct1.Tref_,
220 ct1.Y()/eofs.Y()*ct1.Hsref_
221 + ct2.Y()/eofs.Y()*ct2.Hsref_
222 );
223 }
224}
225
226
227template<class EquationOfState>
228inline Foam::hConstThermo<EquationOfState> Foam::operator*
229(
230 const scalar s,
231 const hConstThermo<EquationOfState>& ct
232)
233{
234 return hConstThermo<EquationOfState>
235 (
236 s*static_cast<const EquationOfState&>(ct),
237 ct.Cp_,
238 ct.Hf_,
239 ct.Tref_,
240 ct.Hsref_
241 );
242}
243
244
245template<class EquationOfState>
246inline Foam::hConstThermo<EquationOfState> Foam::operator==
247(
248 const hConstThermo<EquationOfState>& ct1,
249 const hConstThermo<EquationOfState>& ct2
250)
251{
252 EquationOfState eofs
253 (
254 static_cast<const EquationOfState&>(ct1)
255 == static_cast<const EquationOfState&>(ct2)
256 );
257
258 return hConstThermo<EquationOfState>
259 (
260 eofs,
261 ct2.Y()/eofs.Y()*ct2.Cp_
262 - ct1.Y()/eofs.Y()*ct1.Cp_,
263 ct2.Y()/eofs.Y()*ct2.Hf_
264 - ct1.Y()/eofs.Y()*ct1.Hf_,
265 ct1.Tref_,
266 ct2.Y()/eofs.Y()*ct2.Hsref_
267 - ct1.Y()/eofs.Y()*ct1.Hsref_
268 );
269}
270
271
272// ************************************************************************* //
scalar Hs(const scalar p, const scalar T) const
Definition EtoHthermo.H:17
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Constant properties thermodynamics package templated into the EquationOfState.
scalar limit(const scalar T) const
Limit temperature to be within the range.
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
scalar Hc() const
Chemical enthalpy [J/kg].
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
autoPtr< hConstThermo > clone() const
Construct and return a clone.
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
PtrList< volScalarField > & Y
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))
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
const volScalarField & cp