hRefConstThermoI.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) 2015-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 cp,
33  const scalar hf,
34  const scalar tref,
35  const scalar href
36 )
37 :
38  EquationOfState(st),
39  Cp_(cp),
40  Hf_(hf),
41  Tref_(tref),
42  Href_(href)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class EquationOfState>
50 (
51  const word& name,
52  const hRefConstThermo& ct
53 )
54 :
55  EquationOfState(name, ct),
56  Cp_(ct.Cp_),
57  Hf_(ct.Hf_),
58  Tref_(ct.Tref_),
59  Href_(ct.Href_)
60 {}
61 
62 
63 template<class EquationOfState>
66 {
68  (
70  );
71 }
72 
73 
74 template<class EquationOfState>
77 {
79  (
81  );
82 }
83 
84 
85 template<class EquationOfState>
88 {
90  (
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 template<class EquationOfState>
100 (
101  const scalar T
102 ) const
103 {
104  return T;
105 }
106 
107 
108 template<class EquationOfState>
110 (
111  const scalar p,
112  const scalar T
113 ) const
114 {
115  return Cp_ + EquationOfState::cp(p, T);
116 }
117 
118 
119 template<class EquationOfState>
121 (
122  const scalar p, const scalar T
123 ) const
124 {
125  return Cp_*(T-Tref_) + Href_ + Hf_ + EquationOfState::h(p, T);
126 }
127 
128 
129 template<class EquationOfState>
131 (
132  const scalar p, const scalar T
133 ) const
134 {
135  return Cp_*(T-Tref_) + Href_ + EquationOfState::h(p, T);
136 }
137 
138 
139 template<class EquationOfState>
141 {
142  return Hf_;
143 }
144 
145 
146 template<class EquationOfState>
148 (
149  const scalar p, const scalar T
150 ) const
151 {
152  return Cp_*log(T/Tstd) + EquationOfState::s(p, T);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157 
158 template<class EquationOfState>
159 inline void Foam::hRefConstThermo<EquationOfState>::operator+=
160 (
162 )
163 {
164  scalar molr1 = this->nMoles();
165 
166  EquationOfState::operator+=(ct);
167 
168  molr1 /= this->nMoles();
169  scalar molr2 = ct.nMoles()/this->nMoles();
170 
171  Cp_ = molr1*Cp_ + molr2*ct.Cp_;
172  Hf_ = molr1*Hf_ + molr2*ct.Hf_;
173 }
174 
175 
176 template<class EquationOfState>
177 inline void Foam::hRefConstThermo<EquationOfState>::operator-=
178 (
180 )
181 {
182  scalar molr1 = this->nMoles();
183 
184  EquationOfState::operator-=(ct);
185 
186  molr1 /= this->nMoles();
187  scalar molr2 = ct.nMoles()/this->nMoles();
188 
189  Cp_ = molr1*Cp_ - molr2*ct.Cp_;
190  Hf_ = molr1*Hf_ - molr2*ct.Hf_;
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
195 
196 template<class EquationOfState>
197 inline Foam::hRefConstThermo<EquationOfState> Foam::operator+
198 (
201 )
202 {
203  EquationOfState eofs
204  (
205  static_cast<const EquationOfState&>(ct1)
206  + static_cast<const EquationOfState&>(ct2)
207  );
208 
210  (
211  eofs,
212  ct1.nMoles()/eofs.nMoles()*ct1.Cp_
213  + ct2.nMoles()/eofs.nMoles()*ct2.Cp_,
214  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
215  + ct2.nMoles()/eofs.nMoles()*ct2.Hf_,
216  ct1.nMoles()/eofs.nMoles()*ct1.Tref_
217  + ct2.nMoles()/eofs.nMoles()*ct2.Tref_,
218  ct1.nMoles()/eofs.nMoles()*ct1.Href_
219  + ct2.nMoles()/eofs.nMoles()*ct2.Href_
220  );
221 }
222 
223 
224 template<class EquationOfState>
225 inline Foam::hRefConstThermo<EquationOfState> Foam::operator-
226 (
229 )
230 {
231  EquationOfState eofs
232  (
233  static_cast<const EquationOfState&>(ct1)
234  - static_cast<const EquationOfState&>(ct2)
235  );
236 
238  (
239  eofs,
240  ct1.nMoles()/eofs.nMoles()*ct1.Cp_
241  - ct2.nMoles()/eofs.nMoles()*ct2.Cp_,
242  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
243  - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
244  );
245 }
246 
247 
248 template<class EquationOfState>
249 inline Foam::hRefConstThermo<EquationOfState> Foam::operator*
250 (
251  const scalar s,
253 )
254 {
256  (
257  s*static_cast<const EquationOfState&>(ct),
258  ct.Cp_,
259  ct.Hf_,
260  ct.Tref_,
261  ct.Href_
262  );
263 }
264 
265 
266 template<class EquationOfState>
267 inline Foam::hRefConstThermo<EquationOfState> Foam::operator==
268 (
271 )
272 {
273  return ct2 - ct1;
274 }
275 
276 
277 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
scalar s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
static autoPtr< hRefConstThermo > New(Istream &is)
Selector from Istream.
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].
scalar hc() const
Chemical enthalpy [J/kmol].
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
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))
A class for handling words, derived from string.
Definition: word.H:59
const dimensionedScalar Tstd
Standard temperature.
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
const volScalarField & cp
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
autoPtr< hRefConstThermo > clone() const
Construct and return a clone.
scalar cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [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
Constant properties thermodynamics package templated into the EquationOfState.