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-2016 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 template<class EquationOfState>
82 {
84  (
86  );
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
92 template<class EquationOfState>
94 (
95  const scalar T
96 ) const
97 {
98  return T;
99 }
100 
101 
102 template<class EquationOfState>
103 inline Foam::scalar Foam::eConstThermo<EquationOfState>::cp
104 (
105  const scalar p,
106  const scalar T
107 ) const
108 {
109  return Cv_ + this->cpMcv(p, T) + EquationOfState::cp(p, T);
110 }
111 
112 
113 template<class EquationOfState>
114 inline Foam::scalar Foam::eConstThermo<EquationOfState>::ha
115 (
116  const scalar p,
117  const scalar T
118 ) const
119 {
120  return cp(p, T)*T + Hf_ + EquationOfState::h(p, T);
121 }
122 
123 
124 template<class EquationOfState>
125 inline Foam::scalar Foam::eConstThermo<EquationOfState>::hs
126 (
127  const scalar p,
128  const scalar T
129 ) const
130 {
131  return cp(p, T)*T + EquationOfState::h(p, T);
132 }
133 
134 
135 template<class EquationOfState>
136 inline Foam::scalar Foam::eConstThermo<EquationOfState>::hc() const
137 {
138  return Hf_;
139 }
140 
141 
142 template<class EquationOfState>
143 inline Foam::scalar Foam::eConstThermo<EquationOfState>::s
144 (
145  const scalar p,
146  const scalar T
147 ) const
148 {
149  return cp()*log(T/Tstd) + EquationOfState::s(p, T);
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
155 template<class EquationOfState>
156 inline void Foam::eConstThermo<EquationOfState>::operator+=
157 (
159 )
160 {
161  scalar molr1 = this->nMoles();
162 
163  EquationOfState::operator+=(ct);
164 
165  molr1 /= this->nMoles();
166  scalar molr2 = ct.nMoles()/this->nMoles();
167 
168  Cv_ = molr1*Cv_ + molr2*ct.Cv_;
169  Hf_ = molr1*Hf_ + molr2*ct.Hf_;
170 }
171 
172 
173 template<class EquationOfState>
174 inline void Foam::eConstThermo<EquationOfState>::operator-=
175 (
177 )
178 {
179  scalar molr1 = this->nMoles();
180 
181  EquationOfState::operator-=(ct);
182 
183  molr1 /= this->nMoles();
184  scalar molr2 = ct.nMoles()/this->nMoles();
185 
186  Cv_ = molr1*Cv_ - molr2*ct.Cv_;
187  Hf_ = molr1*Hf_ - molr2*ct.Hf_;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
192 
193 template<class EquationOfState>
194 inline 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 
207  (
208  eofs,
209  ct1.nMoles()/eofs.nMoles()*ct1.Cv_
210  + ct2.nMoles()/eofs.nMoles()*ct2.Cv_,
211  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
212  + ct2.nMoles()/eofs.nMoles()*ct2.Hf_
213  );
214 }
215 
216 
217 template<class EquationOfState>
218 inline Foam::eConstThermo<EquationOfState> Foam::operator-
219 (
222 )
223 {
224  EquationOfState eofs
225  (
226  static_cast<const EquationOfState&>(ct1)
227  - static_cast<const EquationOfState&>(ct2)
228  );
229 
231  (
232  eofs,
233  ct1.nMoles()/eofs.nMoles()*ct1.Cv_
234  - ct2.nMoles()/eofs.nMoles()*ct2.Cv_,
235  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
236  - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
237  );
238 }
239 
240 
241 template<class EquationOfState>
242 inline Foam::eConstThermo<EquationOfState> Foam::operator*
243 (
244  const scalar s,
246 )
247 {
249  (
250  s*static_cast<const EquationOfState&>(ct),
251  ct.Cv_,
252  ct.Hf_
253  );
254 }
255 
256 
257 template<class EquationOfState>
258 inline Foam::eConstThermo<EquationOfState> Foam::operator==
259 (
262 )
263 {
264  return ct2 - ct1;
265 }
266 
267 
268 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar hs(const scalar p, const scalar T) const
Sensible Enthalpy [J/kmol].
autoPtr< eConstThermo > clone() const
Construct and return a clone.
Definition: eConstThermoI.H:59
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/(kmol K)].
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/kmol].
const dimensionedScalar Tstd
Standard temperature.
const volScalarField & cp
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
scalar hc() const
Chemical enthalpy [J/kmol].
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: eConstThermoI.H:94
static autoPtr< eConstThermo > New(Istream &is)
Definition: eConstThermoI.H:70