eRefConstThermoI.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) 2018 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  const scalar tref,
37  const scalar eref
38 )
39 :
40  EquationOfState(st),
41  Cv_(cv),
42  Hf_(hf),
43  Tref_(tref),
44  Eref_(eref)
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class EquationOfState>
52 (
53  const word& name,
54  const eRefConstThermo& ct
55 )
56 :
57  EquationOfState(name, ct),
58  Cv_(ct.Cv_),
59  Hf_(ct.Hf_),
60  Tref_(ct.Tref_),
61  Eref_(ct.Eref_)
62 {}
63 
64 
65 template<class EquationOfState>
68 {
70  (
72  );
73 }
74 
75 
76 template<class EquationOfState>
79 {
81  (
83  );
84 }
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
89 template<class EquationOfState>
91 (
92  const scalar T
93 ) const
94 {
95  return T;
96 }
97 
98 
99 template<class EquationOfState>
101 (
102  const scalar p,
103  const scalar T
104 ) const
105 {
106  return Cv_ + EquationOfState::Cv(p, T);
107 }
108 
109 
110 template<class EquationOfState>
112 (
113  const scalar p, const scalar T
114 ) const
115 {
116  return Cv_*(T - Tref_) + Eref_ + EquationOfState::E(p, T);
117 }
118 
119 
120 template<class EquationOfState>
122 {
123  return Hf_;
124 }
125 
126 
127 template<class EquationOfState>
129 (
130  const scalar p, const scalar T
131 ) const
132 {
133  return Es(p, T) + Hc();
134 }
135 
136 
137 template<class EquationOfState>
139 (
140  const scalar p, const scalar T
141 ) const
142 {
143  return Cp(p, T)*log(T/Tstd) + EquationOfState::S(p, T);
144 }
145 
146 
147 template<class EquationOfState>
149 (
150  const scalar p, const scalar T
151 ) const
152 {
153  return 0;
154 }
155 
156 
157 template<class EquationOfState>
159 (
160  const scalar p, const scalar T
161 ) const
162 {
163  return 0;
164 }
165 
166 
167 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
168 
169 template<class EquationOfState>
171 (
173 )
174 {
175  scalar Y1 = this->Y();
176 
177  EquationOfState::operator+=(ct);
178 
179  if (mag(this->Y()) > SMALL)
180  {
181  Y1 /= this->Y();
182  const scalar Y2 = ct.Y()/this->Y();
183 
184  Cv_ = Y1*Cv_ + Y2*ct.Cv_;
185  Hf_ = Y1*Hf_ + Y2*ct.Hf_;
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
191 
192 template<class EquationOfState>
193 inline Foam::eRefConstThermo<EquationOfState> Foam::operator+
194 (
197 )
198 {
199  EquationOfState eofs
200  (
201  static_cast<const EquationOfState&>(ct1)
202  + static_cast<const EquationOfState&>(ct2)
203  );
204 
205  if (mag(eofs.Y()) < SMALL)
206  {
208  (
209  eofs,
210  ct1.Cv_,
211  ct1.Hf_,
212  ct1.Tref_,
213  ct1.Eref_
214  );
215  }
216  else
217  {
218  return eRefConstThermo<EquationOfState>
219  (
220  eofs,
221  ct1.Y()/eofs.Y()*ct1.Cv_
222  + ct2.Y()/eofs.Y()*ct2.Cv_,
223  ct1.Y()/eofs.Y()*ct1.Hf_
224  + ct2.Y()/eofs.Y()*ct2.Hf_,
225  ct1.Y()/eofs.Y()*ct1.Tref_
226  + ct2.Y()/eofs.Y()*ct2.Tref_,
227  ct1.Y()/eofs.Y()*ct1.Eref_
228  + ct2.Y()/eofs.Y()*ct2.Eref_
229  );
230  }
231 }
232 
233 
234 template<class EquationOfState>
235 inline Foam::eRefConstThermo<EquationOfState> Foam::operator*
236 (
237  const scalar s,
238  const eRefConstThermo<EquationOfState>& ct
239 )
240 {
241  return eRefConstThermo<EquationOfState>
242  (
243  s*static_cast<const EquationOfState&>(ct),
244  ct.Cv_,
245  ct.Hf_,
246  ct.Tref_,
247  ct.Eref_
248  );
249 }
250 
251 
252 template<class EquationOfState>
253 inline Foam::eRefConstThermo<EquationOfState> Foam::operator==
254 (
255  const eRefConstThermo<EquationOfState>& ct1,
256  const eRefConstThermo<EquationOfState>& ct2
257 )
258 {
259  EquationOfState eofs
260  (
261  static_cast<const EquationOfState&>(ct1)
262  == static_cast<const EquationOfState&>(ct2)
263  );
264 
265  return eRefConstThermo<EquationOfState>
266  (
267  eofs,
268  ct2.Y()/eofs.Y()*ct2.Cv_
269  - ct1.Y()/eofs.Y()*ct1.Cv_,
270  ct2.Y()/eofs.Y()*ct2.Hf_
271  - ct1.Y()/eofs.Y()*ct1.Hf_
272  );
273 }
274 
275 
276 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::eRefConstThermo::dCpdT
scalar dCpdT(const scalar p, const scalar T) const
Temperature derivative of heat capacity at constant pressure.
Definition: eRefConstThermoI.H:159
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
Cv
const volScalarField & Cv
Definition: EEqn.H:8
Foam::eRefConstThermo::Cv
scalar Cv(const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
Definition: eRefConstThermoI.H:101
Foam::eRefConstThermo::limit
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
Definition: eRefConstThermoI.H:91
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
Foam::eRefConstThermo::dGdT
scalar dGdT(const scalar p, const scalar T) const
Derivative of Gibbs free energy w.r.t. temperature.
Definition: eRefConstThermoI.H:149
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::eRefConstThermo::Es
scalar Es(const scalar p, const scalar T) const
Sensible internal energy [J/kg].
Definition: eRefConstThermoI.H:112
Foam::eRefConstThermo::clone
autoPtr< eRefConstThermo > clone() const
Construct and return a clone.
Definition: eRefConstThermoI.H:67
Foam::eRefConstThermo::Hc
scalar Hc() const
Chemical enthalpy [J/kg].
Definition: eRefConstThermoI.H:121
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::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::eRefConstThermo::S
scalar S(const scalar p, const scalar T) const
Entropy [J/(kg K)].
Definition: eRefConstThermoI.H:139
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::eRefConstThermo
Constant properties thermodynamics package templated into the EquationOfState.
Definition: eRefConstThermo.H:48
Foam::eRefConstThermo::Ea
scalar Ea(const scalar p, const scalar T) const
Absolute internal energy [J/kg].
Definition: eRefConstThermoI.H:129
Foam::eRefConstThermo::New
static autoPtr< eRefConstThermo > New(const dictionary &dict)
Selector from dictionary.
Definition: eRefConstThermoI.H:78