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 | 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 cv,
35  const scalar hf
36 )
37 :
38  EquationOfState(st),
39  Cv_(cv),
40  Hf_(hf)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class EquationOfState>
48 (
49  const word& name,
50  const eConstThermo<EquationOfState>& ct
51 )
52 :
53  EquationOfState(name, ct),
54  Cv_(ct.Cv_),
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 Cv_ + EquationOfState::Cv(p, T);
95 }
96 
97 
98 template<class EquationOfState>
100 (
101  const scalar p,
102  const scalar T
103 ) const
104 {
105  return Cv_*T + EquationOfState::E(p, T);
106 }
107 
108 
109 template<class EquationOfState>
110 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Ea
111 (
112  const scalar p,
113  const scalar T
114 ) const
115 {
116  return Es(p, T) + Hc();
117 }
118 
119 
120 template<class EquationOfState>
121 inline Foam::scalar Foam::eConstThermo<EquationOfState>::Hc() const
122 {
123  return Hf_;
124 }
125 
126 
127 template<class EquationOfState>
128 inline Foam::scalar Foam::eConstThermo<EquationOfState>::S
129 (
130  const scalar p,
131  const scalar T
132 ) const
133 {
134  return Cp(p, T)*log(T/Tstd) + EquationOfState::S(p, T);
135 }
136 
137 
138 template<class EquationOfState>
140 (
141  const scalar p,
142  const scalar T
143 ) const
144 {
145  return 0;
146 }
147 
148 
149 template<class EquationOfState>
151 (
152  const scalar p,
153  const scalar T
154 ) const
155 {
156  return 0;
157 }
158 
159 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
160 
161 template<class EquationOfState>
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  const scalar Y2 = ct.Y()/this->Y();
175 
176  Cv_ = Y1*Cv_ + Y2*ct.Cv_;
177  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
178  }
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
183 
184 template<class EquationOfState>
185 inline Foam::eConstThermo<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.Cv_,
203  ct1.Hf_
204  );
205  }
206  else
207  {
208  return eConstThermo<EquationOfState>
209  (
210  eofs,
211  ct1.Y()/eofs.Y()*ct1.Cv_
212  + ct2.Y()/eofs.Y()*ct2.Cv_,
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::eConstThermo<EquationOfState> Foam::operator*
222 (
223  const scalar s,
224  const eConstThermo<EquationOfState>& ct
225 )
226 {
227  return eConstThermo<EquationOfState>
228  (
229  s*static_cast<const EquationOfState&>(ct),
230  ct.Cv_,
231  ct.Hf_
232  );
233 }
234 
235 
236 template<class EquationOfState>
237 inline Foam::eConstThermo<EquationOfState> Foam::operator==
238 (
239  const eConstThermo<EquationOfState>& ct1,
240  const eConstThermo<EquationOfState>& ct2
241 )
242 {
243  EquationOfState eofs
244  (
245  static_cast<const EquationOfState&>(ct1)
246  == static_cast<const EquationOfState&>(ct2)
247  );
248 
249  return eConstThermo<EquationOfState>
250  (
251  eofs,
252  ct2.Y()/eofs.Y()*ct2.Cv_
253  - ct1.Y()/eofs.Y()*ct1.Cv_,
254  ct2.Y()/eofs.Y()*ct2.Hf_
255  - ct1.Y()/eofs.Y()*ct1.Hf_
256  );
257 }
258 
259 
260 // ************************************************************************* //
Foam::eConstThermo::Cv
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
Definition: eConstThermoI.H:89
Foam::eConstThermo::clone
autoPtr< eConstThermo > clone() const
Construct and return a clone.
Definition: eConstThermoI.H:61
Foam::eConstThermo::New
static autoPtr< eConstThermo > New(const dictionary &dict)
Definition: eConstThermoI.H:69
p
volScalarField & p
Definition: createFieldRefs.H:8
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
Foam::eConstThermo
Constant properties thermodynamics package templated on an equation of state.
Definition: eConstThermo.H:53
Foam::eConstThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: eConstThermoI.H:79
Foam::eConstThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: eConstThermoI.H:140
Cv
const volScalarField & Cv
Definition: EEqn.H:8
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
Foam::eConstThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: eConstThermoI.H:121
Foam::eConstThermo::Ea
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
Definition: eConstThermoI.H:111
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
Foam::eConstThermo::Es
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
Definition: eConstThermoI.H:100
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::eConstThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: eConstThermoI.H:129
Es
scalar Es(const scalar p, const scalar T) const
Definition: HtoEthermo.H:17
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::eConstThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: eConstThermoI.H:151