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 | 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 cp,
33  const scalar hf
34 )
35 :
36  EquationOfState(st),
37  Cp_(cp),
38  Hf_(hf)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class EquationOfState>
46 (
47  const word& name,
48  const hConstThermo& ct
49 )
50 :
51  EquationOfState(name, ct),
52  Cp_(ct.Cp_),
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::hConstThermo<EquationOfState>::cp
104 (
105  const scalar p,
106  const scalar T
107 ) const
108 {
109  return Cp_ + EquationOfState::cp(p, T);
110 }
111 
112 
113 template<class EquationOfState>
114 inline Foam::scalar Foam::hConstThermo<EquationOfState>::ha
115 (
116  const scalar p, const scalar T
117 ) const
118 {
119  return Cp_*T + Hf_ + EquationOfState::h(p, T);
120 }
121 
122 
123 template<class EquationOfState>
124 inline Foam::scalar Foam::hConstThermo<EquationOfState>::hs
125 (
126  const scalar p, const scalar T
127 ) const
128 {
129  return Cp_*T + EquationOfState::h(p, T);
130 }
131 
132 
133 template<class EquationOfState>
134 inline Foam::scalar Foam::hConstThermo<EquationOfState>::hc() const
135 {
136  return Hf_;
137 }
138 
139 
140 template<class EquationOfState>
141 inline Foam::scalar Foam::hConstThermo<EquationOfState>::s
142 (
143  const scalar p, const scalar T
144 ) const
145 {
146  return Cp_*log(T/Tstd) + EquationOfState::s(p, T);
147 }
148 
149 
150 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
151 
152 template<class EquationOfState>
153 inline void Foam::hConstThermo<EquationOfState>::operator+=
154 (
156 )
157 {
158  scalar molr1 = this->nMoles();
159 
160  EquationOfState::operator+=(ct);
161 
162  molr1 /= this->nMoles();
163  scalar molr2 = ct.nMoles()/this->nMoles();
164 
165  Cp_ = molr1*Cp_ + molr2*ct.Cp_;
166  Hf_ = molr1*Hf_ + molr2*ct.Hf_;
167 }
168 
169 
170 template<class EquationOfState>
171 inline void Foam::hConstThermo<EquationOfState>::operator-=
172 (
174 )
175 {
176  scalar molr1 = this->nMoles();
177 
178  EquationOfState::operator-=(ct);
179 
180  molr1 /= this->nMoles();
181  scalar molr2 = ct.nMoles()/this->nMoles();
182 
183  Cp_ = molr1*Cp_ - molr2*ct.Cp_;
184  Hf_ = molr1*Hf_ - molr2*ct.Hf_;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
189 
190 template<class EquationOfState>
191 inline Foam::hConstThermo<EquationOfState> Foam::operator+
192 (
195 )
196 {
197  EquationOfState eofs
198  (
199  static_cast<const EquationOfState&>(ct1)
200  + static_cast<const EquationOfState&>(ct2)
201  );
202 
204  (
205  eofs,
206  ct1.nMoles()/eofs.nMoles()*ct1.Cp_
207  + ct2.nMoles()/eofs.nMoles()*ct2.Cp_,
208  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
209  + ct2.nMoles()/eofs.nMoles()*ct2.Hf_
210  );
211 }
212 
213 
214 template<class EquationOfState>
215 inline Foam::hConstThermo<EquationOfState> Foam::operator-
216 (
219 )
220 {
221  EquationOfState eofs
222  (
223  static_cast<const EquationOfState&>(ct1)
224  - static_cast<const EquationOfState&>(ct2)
225  );
226 
228  (
229  eofs,
230  ct1.nMoles()/eofs.nMoles()*ct1.Cp_
231  - ct2.nMoles()/eofs.nMoles()*ct2.Cp_,
232  ct1.nMoles()/eofs.nMoles()*ct1.Hf_
233  - ct2.nMoles()/eofs.nMoles()*ct2.Hf_
234  );
235 }
236 
237 
238 template<class EquationOfState>
239 inline Foam::hConstThermo<EquationOfState> Foam::operator*
240 (
241  const scalar s,
243 )
244 {
246  (
247  s*static_cast<const EquationOfState&>(ct),
248  ct.Cp_,
249  ct.Hf_
250  );
251 }
252 
253 
254 template<class EquationOfState>
255 inline Foam::hConstThermo<EquationOfState> Foam::operator==
256 (
259 )
260 {
261  return ct2 - ct1;
262 }
263 
264 
265 // ************************************************************************* //
dictionary dict
Constant properties thermodynamics package templated into the EquationOfState.
Definition: hConstThermo.H:46
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].
static autoPtr< hConstThermo > New(Istream &is)
Selector from Istream.
Definition: hConstThermoI.H:70
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 cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kmol K)].
A class for handling words, derived from string.
Definition: word.H:59
const dimensionedScalar Tstd
Standard temperature.
scalar ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kmol].
const volScalarField & cp
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
autoPtr< hConstThermo > clone() const
Construct and return a clone.
Definition: hConstThermoI.H:59
scalar hc() const
Chemical enthalpy [J/kmol].
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 s(const scalar p, const scalar T) const
Entropy [J/(kmol K)].
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hConstThermoI.H:94