Loading...
Searching...
No Matches
eConstThermoI.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::eConstThermo<EquationOfState>::eConstThermo
32(
33 const EquationOfState& st,
34 const scalar cv,
35 const scalar hf,
36 const scalar Tref,
37 const scalar Esref
38)
39:
40 EquationOfState(st),
41 Cv_(cv),
42 Hf_(hf),
43 Tref_(Tref),
44 Esref_(Esref)
45{}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
50template<class EquationOfState>
51inline Foam::eConstThermo<EquationOfState>::eConstThermo
52(
53 const word& name,
55)
56:
57 EquationOfState(name, ct),
58 Cv_(ct.Cv_),
59 Hf_(ct.Hf_),
60 Tref_(ct.Tref_),
61 Esref_(ct.Esref_)
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 Cv_ + EquationOfState::Cv(p, T);
101}
102
103
104template<class EquationOfState>
106(
107 const scalar p,
108 const scalar T
109) const
110{
111 return Cv_*(T - Tref_) + Esref_ + EquationOfState::E(p, T);
112}
113
114
115template<class EquationOfState>
117(
118 const scalar p,
119 const scalar T
120) const
121{
122 return Es(p, T) + Hc();
123}
124
125
126template<class EquationOfState>
128{
129 return Hf_;
130}
131
132
133template<class EquationOfState>
135(
136 const scalar p,
137 const scalar T
138) const
139{
140 return Cp(p, T)*log(T/Tstd) + EquationOfState::S(p, T);
141}
142
143
144template<class EquationOfState>
146(
147 const scalar T
148) const
149{
150 return
151 Cv_*(T - Tref_) + Esref_ + Hc() + Pstd/EquationOfState::rho(Pstd, T)
152 - S(Pstd, T)*T;
153}
154
155
156template<class EquationOfState>
158(
159 const scalar p,
160 const scalar T
161) const
162{
164 return 0;
165}
166
167// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
168
169template<class EquationOfState>
171(
173)
174{
175 scalar Y1 = this->Y();
176
177 EquationOfState::operator+=(ct);
178
179 if (mag(this->Y()) > SMALL)
180 {
181 Y1 /= this->Y();
182 const scalar Y2 = ct.Y()/this->Y();
183
184 Cv_ = Y1*Cv_ + Y2*ct.Cv_;
185 Hf_ = Y1*Hf_ + Y2*ct.Hf_;
186 Esref_ = Y1*Esref_ + Y2*ct.Esref_;
187 }
188}
189
190
191// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
192
193template<class EquationOfState>
194inline Foam::eConstThermo<EquationOfState> Foam::operator+
195(
198)
199{
200 EquationOfState eofs
201 (
202 static_cast<const EquationOfState&>(ct1)
203 + static_cast<const EquationOfState&>(ct2)
204 );
205
206 if (mag(eofs.Y()) < SMALL)
207 {
209 (
210 eofs,
211 ct1.Cv_,
212 ct1.Hf_,
213 ct1.Tref_,
214 ct1.Esref_
215 );
216 }
217 else
218 {
220 (
221 eofs,
222 ct1.Y()/eofs.Y()*ct1.Cv_
223 + ct2.Y()/eofs.Y()*ct2.Cv_,
224 ct1.Y()/eofs.Y()*ct1.Hf_
225 + ct2.Y()/eofs.Y()*ct2.Hf_,
226 ct1.Tref_,
227 ct1.Y()/eofs.Y()*ct1.Esref_
228 + ct2.Y()/eofs.Y()*ct2.Esref_
229 );
230 }
231}
232
233
234template<class EquationOfState>
235inline Foam::eConstThermo<EquationOfState> Foam::operator*
236(
237 const scalar s,
239)
240{
242 (
243 s*static_cast<const EquationOfState&>(ct),
244 ct.Cv_,
245 ct.Hf_,
246 ct.Tref_,
247 ct.Esref_
248 );
249}
250
251
252template<class EquationOfState>
253inline Foam::eConstThermo<EquationOfState> Foam::operator==
254(
257)
258{
259 EquationOfState eofs
260 (
261 static_cast<const EquationOfState&>(ct1)
262 == static_cast<const EquationOfState&>(ct2)
263 );
264
266 (
267 eofs,
268 ct2.Y()/eofs.Y()*ct2.Cv_
269 - ct1.Y()/eofs.Y()*ct1.Cv_,
270 ct2.Y()/eofs.Y()*ct2.Hf_
271 - ct1.Y()/eofs.Y()*ct1.Hf_,
272 ct1.Tref_,
273 ct2.Y()/eofs.Y()*ct2.Esref_
274 - ct1.Y()/eofs.Y()*ct1.Esref_
275 );
276}
277
278
279// ************************************************************************* //
scalar Es(const scalar p, const scalar T) const
Definition HtoEthermo.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 on an equation of state.
scalar limit(const scalar T) const
Limit temperature to be within the range.
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
scalar Hc() const
Chemical enthalpy [J/kg].
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
scalar Gstd(const scalar T) const
Gibbs free energy of the mixture in the standard state [J/kg].
static autoPtr< eConstThermo > New(const dictionary &dict)
autoPtr< eConstThermo > clone() const
Construct and return a clone.
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & Cp
Definition EEqn.H:7
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
auto & name
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))
const scalar Pstd
Standard pressure: default in [Pa].
const scalar Tstd
Standard temperature: default in [K].
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