hConstThermoI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 // * * * * * * * * * * * * * * * 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 Cp_ + EquationOfState::Cp(p, T);
99 }
100 
101 
102 template<class EquationOfState>
103 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Ha
104 (
105  const scalar p, const scalar T
106 ) const
107 {
108  return Cp_*T + Hf_ + EquationOfState::H(p, T);
109 }
110 
111 
112 template<class EquationOfState>
113 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hs
114 (
115  const scalar p, const scalar T
116 ) const
117 {
118  return Cp_*T + EquationOfState::H(p, T);
119 }
120 
121 
122 template<class EquationOfState>
123 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hc() const
124 {
125  return Hf_;
126 }
127 
128 
129 template<class EquationOfState>
130 inline Foam::scalar Foam::hConstThermo<EquationOfState>::S
131 (
132  const scalar p, const scalar T
133 ) const
134 {
135  return Cp_*log(T/Tstd) + EquationOfState::S(p, T);
136 }
137 
138 
139 template<class EquationOfState>
141 (
142  const scalar p, const scalar T
143 ) const
144 {
145  return 0;
146 }
147 
148 
149 template<class EquationOfState>
151 (
152  const scalar p, const scalar T
153 ) const
154 {
155  return 0;
156 }
157 
158 
159 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
160 
161 template<class EquationOfState>
162 inline void Foam::hConstThermo<EquationOfState>::operator+=
163 (
165 )
166 {
167  scalar Y1 = this->Y();
168 
169  EquationOfState::operator+=(ct);
170 
171  if (mag(this->Y()) > small)
172  {
173  Y1 /= this->Y();
174  scalar Y2 = ct.Y()/this->Y();
175 
176  Cp_ = Y1*Cp_ + Y2*ct.Cp_;
177  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
178  }
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
183 
184 template<class EquationOfState>
185 inline Foam::hConstThermo<EquationOfState> Foam::operator+
186 (
189 )
190 {
191  EquationOfState eofs
192  (
193  static_cast<const EquationOfState&>(ct1)
194  + static_cast<const EquationOfState&>(ct2)
195  );
196 
197  if (mag(eofs.Y()) < small)
198  {
200  (
201  eofs,
202  ct1.Cp_,
203  ct1.Hf_
204  );
205  }
206  else
207  {
209  (
210  eofs,
211  ct1.Y()/eofs.Y()*ct1.Cp_
212  + ct2.Y()/eofs.Y()*ct2.Cp_,
213  ct1.Y()/eofs.Y()*ct1.Hf_
214  + ct2.Y()/eofs.Y()*ct2.Hf_
215  );
216  }
217 }
218 
219 
220 template<class EquationOfState>
221 inline Foam::hConstThermo<EquationOfState> Foam::operator*
222 (
223  const scalar s,
225 )
226 {
228  (
229  s*static_cast<const EquationOfState&>(ct),
230  ct.Cp_,
231  ct.Hf_
232  );
233 }
234 
235 
236 template<class EquationOfState>
237 inline Foam::hConstThermo<EquationOfState> Foam::operator==
238 (
241 )
242 {
243  EquationOfState eofs
244  (
245  static_cast<const EquationOfState&>(ct1)
246  == static_cast<const EquationOfState&>(ct2)
247  );
248 
250  (
251  eofs,
252  ct2.Y()/eofs.Y()*ct2.Cp_
253  - ct1.Y()/eofs.Y()*ct1.Cp_,
254  ct2.Y()/eofs.Y()*ct2.Hf_
255  - ct1.Y()/eofs.Y()*ct1.Hf_
256  );
257 }
258 
259 
260 // ************************************************************************* //
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
dictionary dict
Constant properties thermodynamics package templated into the EquationOfState.
Definition: hConstThermo.H:46
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hConstThermoI.H:83
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
autoPtr< hConstThermo > clone() const
Construct and return a clone.
Definition: hConstThermoI.H:59
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: hConstThermoI.H:70
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
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
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: hConstThermoI.H:93
const dimensionedScalar Tstd
Standard temperature.
const volScalarField & cp
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
scalar Hc() const
Chemical enthalpy [J/kg].