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 | Copyright (C) 2011-2017 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 template<class EquationOfState>
30 (
31  const EquationOfState& st,
32  const scalar cv,
33  const scalar hf
34 )
35 :
36  EquationOfState(st),
37  Cv_(cv),
38  Hf_(hf)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class EquationOfState>
46 (
47  const word& name,
48  const eConstThermo<EquationOfState>& ct
49 )
50 :
51  EquationOfState(name, ct),
52  Cv_(ct.Cv_),
53  Hf_(ct.Hf_)
54 {}
55 
56 
57 template<class EquationOfState>
60 {
62  (
64  );
65 }
66 
67 
68 template<class EquationOfState>
71 {
73  (
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
81 template<class EquationOfState>
83 (
84  const scalar T
85 ) const
86 {
87  return T;
88 }
89 
90 
91 template<class EquationOfState>
93 (
94  const scalar p,
95  const scalar T
96 ) const
97 {
98  return Cv_ + this->CpMCv(p, T) + EquationOfState::Cp(p, T);
99 }
100 
101 
102 template<class EquationOfState>
103 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Ha
104 (
105  const scalar p,
106  const scalar T
107 ) const
108 {
109  return Cp(p, T)*T + Hf_ + EquationOfState::H(p, T);
110 }
111 
112 
113 template<class EquationOfState>
114 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Hs
115 (
116  const scalar p,
117  const scalar T
118 ) const
119 {
120  return Cp(p, T)*T + EquationOfState::H(p, T);
121 }
122 
123 
124 template<class EquationOfState>
125 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Hc() const
126 {
127  return Hf_;
128 }
129 
130 
131 template<class EquationOfState>
132 inline Foam::scalar Foam::eConstThermo<EquationOfState>::S
133 (
134  const scalar p,
135  const scalar T
136 ) const
137 {
138  return Cp(p, T)*log(T/Tstd) + EquationOfState::S(p, T);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
143 
144 template<class EquationOfState>
145 inline void Foam::eConstThermo<EquationOfState>::operator+=
146 (
148 )
149 {
150  scalar Y1 = this->Y();
151 
152  EquationOfState::operator+=(ct);
153 
154  if (mag(this->Y()) > SMALL)
155  {
156  Y1 /= this->Y();
157  const scalar Y2 = ct.Y()/this->Y();
158 
159  Cv_ = Y1*Cv_ + Y2*ct.Cv_;
160  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
161  }
162 }
163 
164 
165 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
166 
167 template<class EquationOfState>
168 inline Foam::eConstThermo<EquationOfState> Foam::operator+
169 (
172 )
173 {
174  EquationOfState eofs
175  (
176  static_cast<const EquationOfState&>(ct1)
177  + static_cast<const EquationOfState&>(ct2)
178  );
179 
180  if (mag(eofs.Y()) < SMALL)
181  {
183  (
184  eofs,
185  ct1.Cv_,
186  ct1.Hf_
187  );
188  }
189  else
190  {
192  (
193  eofs,
194  ct1.Y()/eofs.Y()*ct1.Cv_
195  + ct2.Y()/eofs.Y()*ct2.Cv_,
196  ct1.Y()/eofs.Y()*ct1.Hf_
197  + ct2.Y()/eofs.Y()*ct2.Hf_
198  );
199  }
200 }
201 
202 
203 template<class EquationOfState>
204 inline Foam::eConstThermo<EquationOfState> Foam::operator*
205 (
206  const scalar s,
208 )
209 {
211  (
212  s*static_cast<const EquationOfState&>(ct),
213  ct.Cv_,
214  ct.Hf_
215  );
216 }
217 
218 
219 template<class EquationOfState>
220 inline Foam::eConstThermo<EquationOfState> Foam::operator==
221 (
224 )
225 {
226  EquationOfState eofs
227  (
228  static_cast<const EquationOfState&>(ct1)
229  == static_cast<const EquationOfState&>(ct2)
230  );
231 
233  (
234  eofs,
235  ct2.Y()/eofs.Y()*ct2.Cv_
236  - ct1.Y()/eofs.Y()*ct1.Cv_,
237  ct2.Y()/eofs.Y()*ct2.Hf_
238  - ct1.Y()/eofs.Y()*ct1.Hf_
239  );
240 }
241 
242 
243 // ************************************************************************* //
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: eConstThermoI.H:83
PtrList< volScalarField > & Y1
Definition: YEqns.H:8
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
scalar Hc() const
Chemical enthalpy [J/kg].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
PtrList< volScalarField > & Y2
Definition: YEqns.H:9
Constant properties thermodynamics package templated on an equation of state.
Definition: eConstThermo.H:48
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: eConstThermoI.H:93
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 Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
const dimensionedScalar Tstd
Standard temperature.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< eConstThermo > clone() const
Construct and return a clone.
Definition: eConstThermoI.H:59
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static autoPtr< eConstThermo > New(const dictionary &dict)
Definition: eConstThermoI.H:70
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].