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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class EquationOfState>
32 (
33  const EquationOfState& st,
34  const scalar cp,
35  const scalar hf
36 )
37 :
38  EquationOfState(st),
39  Cp_(cp),
40  Hf_(hf)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class EquationOfState>
48 (
49  const word& name,
50  const hConstThermo& ct
51 )
52 :
53  EquationOfState(name, ct),
54  Cp_(ct.Cp_),
55  Hf_(ct.Hf_)
56 {}
57 
58 
59 template<class EquationOfState>
62 {
64 }
65 
66 
67 template<class EquationOfState>
70 {
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
77 template<class EquationOfState>
79 (
80  const scalar T
81 ) const
82 {
83  return T;
84 }
85 
86 
87 template<class EquationOfState>
89 (
90  const scalar p,
91  const scalar T
92 ) const
93 {
94  return Cp_ + EquationOfState::Cp(p, T);
95 }
96 
97 
98 template<class EquationOfState>
100 (
101  const scalar p, const scalar T
102 ) const
103 {
104  return Hs(p, T) + Hc();
105 }
106 
107 
108 template<class EquationOfState>
109 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hs
110 (
111  const scalar p, const scalar T
112 ) const
113 {
114  return Cp_*T + EquationOfState::H(p, T);
115 }
116 
117 
118 template<class EquationOfState>
119 inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hc() const
120 {
121  return Hf_;
122 }
123 
124 
125 template<class EquationOfState>
126 inline Foam::scalar Foam::hConstThermo<EquationOfState>::S
127 (
128  const scalar p, const scalar T
129 ) const
130 {
131  return Cp_*log(T/Tstd) + EquationOfState::S(p, T);
132 }
133 
134 
135 template<class EquationOfState>
137 (
138  const scalar p, const scalar T
139 ) const
140 {
141  return 0;
142 }
143 
144 
145 template<class EquationOfState>
147 (
148  const scalar p, const scalar T
149 ) const
150 {
151  return 0;
152 }
153 
154 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
155 
156 template<class EquationOfState>
158 (
160 )
161 {
162  scalar Y1 = this->Y();
163 
164  EquationOfState::operator+=(ct);
165 
166  if (mag(this->Y()) > SMALL)
167  {
168  Y1 /= this->Y();
169  scalar Y2 = ct.Y()/this->Y();
170 
171  Cp_ = Y1*Cp_ + Y2*ct.Cp_;
172  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
173  }
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
178 
179 template<class EquationOfState>
180 inline Foam::hConstThermo<EquationOfState> Foam::operator+
181 (
184 )
185 {
186  EquationOfState eofs
187  (
188  static_cast<const EquationOfState&>(ct1)
189  + static_cast<const EquationOfState&>(ct2)
190  );
191 
192  if (mag(eofs.Y()) < SMALL)
193  {
195  (
196  eofs,
197  ct1.Cp_,
198  ct1.Hf_
199  );
200  }
201  else
202  {
203  return hConstThermo<EquationOfState>
204  (
205  eofs,
206  ct1.Y()/eofs.Y()*ct1.Cp_
207  + ct2.Y()/eofs.Y()*ct2.Cp_,
208  ct1.Y()/eofs.Y()*ct1.Hf_
209  + ct2.Y()/eofs.Y()*ct2.Hf_
210  );
211  }
212 }
213 
214 
215 template<class EquationOfState>
216 inline Foam::hConstThermo<EquationOfState> Foam::operator*
217 (
218  const scalar s,
219  const hConstThermo<EquationOfState>& ct
220 )
221 {
222  return hConstThermo<EquationOfState>
223  (
224  s*static_cast<const EquationOfState&>(ct),
225  ct.Cp_,
226  ct.Hf_
227  );
228 }
229 
230 
231 template<class EquationOfState>
232 inline Foam::hConstThermo<EquationOfState> Foam::operator==
233 (
234  const hConstThermo<EquationOfState>& ct1,
235  const hConstThermo<EquationOfState>& ct2
236 )
237 {
238  EquationOfState eofs
239  (
240  static_cast<const EquationOfState&>(ct1)
241  == static_cast<const EquationOfState&>(ct2)
242  );
243 
244  return hConstThermo<EquationOfState>
245  (
246  eofs,
247  ct2.Y()/eofs.Y()*ct2.Cp_
248  - ct1.Y()/eofs.Y()*ct1.Cp_,
249  ct2.Y()/eofs.Y()*ct2.Hf_
250  - ct1.Y()/eofs.Y()*ct1.Hf_
251  );
252 }
253 
254 
255 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::hConstThermo::Hs
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Definition: hConstThermoI.H:110
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
cp
const volScalarField & cp
Definition: setRegionSolidFields.H:8
Hs
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:17
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::hConstThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: hConstThermoI.H:119
Foam::hConstThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: hConstThermoI.H:79
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
Foam::hConstThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: hConstThermoI.H:137
Foam::hConstThermo::New
static autoPtr< hConstThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: hConstThermoI.H:69
Foam::hConstThermo::clone
autoPtr< hConstThermo > clone() const
Construct and return a clone.
Definition: hConstThermoI.H:61
Foam::hConstThermo::Ha
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
Definition: hConstThermoI.H:100
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::hConstThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: hConstThermoI.H:147
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::hConstThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: hConstThermoI.H:127
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::hConstThermo::Cp
scalar Cp(const scalar p, const scalar T) const
Heat capacity at constant pressure [J/(kg K)].
Definition: hConstThermoI.H:89
Foam::hConstThermo
Constant properties thermodynamics package templated into the EquationOfState.
Definition: hConstThermo.H:51