liquidPropertiesI.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 inline Foam::scalar Foam::liquidProperties::limit(const scalar T) const
27 {
28  return T;
29 }
30 
31 
32 inline Foam::scalar Foam::liquidProperties::Y() const
33 {
34  return 1;
35 }
36 
37 
38 inline Foam::scalar Foam::liquidProperties::Tc() const
39 {
40  return Tc_;
41 }
42 
43 
44 inline Foam::scalar Foam::liquidProperties::Pc() const
45 {
46  return Pc_;
47 }
48 
49 
50 inline Foam::scalar Foam::liquidProperties::Vc() const
51 {
52  return Vc_;
53 }
54 
55 
56 inline Foam::scalar Foam::liquidProperties::Zc() const
57 {
58  return Zc_;
59 }
60 
61 
62 inline Foam::scalar Foam::liquidProperties::Tt() const
63 {
64  return Tt_;
65 }
66 
67 
68 inline Foam::scalar Foam::liquidProperties::Pt() const
69 {
70  return Pt_;
71 }
72 
73 
74 inline Foam::scalar Foam::liquidProperties::Tb() const
75 {
76  return Tb_;
77 }
78 
79 
80 inline Foam::scalar Foam::liquidProperties::dipm() const
81 {
82  return dipm_;
83 }
84 
85 
86 inline Foam::scalar Foam::liquidProperties::omega() const
87 {
88  return omega_;
89 }
90 
91 
92 inline Foam::scalar Foam::liquidProperties::delta() const
93 {
94  return delta_;
95 }
96 
97 
98 inline Foam::scalar Foam::liquidProperties::psi(scalar p, scalar T) const
99 {
100  return 0;
101 }
102 
103 
104 inline Foam::scalar Foam::liquidProperties::CpMCv(scalar p, scalar T) const
105 {
106  return 0;
107 }
108 
109 
110 inline Foam::scalar Foam::liquidProperties::Ha(scalar p, scalar T) const
111 {
112  return h(p, T);
113 }
114 
115 
116 inline Foam::scalar Foam::liquidProperties::Hs(scalar p, scalar T) const
117 {
118  return h(p, T);
119 }
120 
121 
122 inline Foam::scalar Foam::liquidProperties::Hc() const
123 {
124  return 0;
125 }
126 
127 
128 inline Foam::scalar Foam::liquidProperties::alphah(scalar p, scalar T) const
129 {
130  return kappa(p, T)/Cp(p, T);
131 }
132 
133 
134 template<class Func>
136 (
137  Func& f,
138  const word& name,
139  const dictionary& dict
140 )
141 {
142  if (dict.found(name))
143  {
144  f = Func(dict.subDict(name));
145  }
146 }
147 
148 
149 template<class Liquid>
151 (
152  Liquid& l,
153  const dictionary& dict
154 )
155 {
156  l.liquidProperties::readIfPresent(dict);
157  readIfPresent(l.rho_, "rho", dict);
158  readIfPresent(l.pv_, "pv", dict);
159  readIfPresent(l.hl_, "hl", dict);
160  readIfPresent(l.Cp_, "Cp", dict);
161  readIfPresent(l.h_, "h", dict);
162  readIfPresent(l.Cpg_, "Cpg", dict);
163  readIfPresent(l.B_, "B", dict);
164  readIfPresent(l.mu_, "mu", dict);
165  readIfPresent(l.mug_, "mug", dict);
166  readIfPresent(l.kappa_, "K", dict);
167  readIfPresent(l.kappag_, "kappag", dict);
168  readIfPresent(l.sigma_, "sigma", dict);
169  readIfPresent(l.D_, "D", dict);
170 }
171 
172 
173 template<class Liquid>
175 (
176  const Liquid& l,
177  Ostream& os
178 ) const
179 {
180  l.liquidProperties::writeData(os); os << nl;
181  l.rho_.writeData(os); os << nl;
182  l.pv_.writeData(os); os << nl;
183  l.hl_.writeData(os); os << nl;
184  l.Cp_.writeData(os); os << nl;
185  l.h_.writeData(os); os << nl;
186  l.Cpg_.writeData(os); os << nl;
187  l.B_.writeData(os); os << nl;
188  l.mu_.writeData(os); os << nl;
189  l.mug_.writeData(os); os << nl;
190  l.kappa_.writeData(os); os << nl;
191  l.kappag_.writeData(os); os << nl;
192  l.sigma_.writeData(os); os << nl;
193  l.D_.writeData(os); os << endl;
194 }
195 
196 
197 // ************************************************************************* //
scalar delta() const
Solubility parameter [(J/m^3)^(1/2)].
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void readIfPresent(const dictionary &dict)
Read and set the properties present it the given dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
scalar Y() const
No of moles of this species in mixture.
scalar Ha(const scalar p, const scalar T) const
Absolute Enthalpy [J/kg].
scalar Pt() const
Triple point pressure [Pa].
virtual scalar kappa(scalar p, scalar T) const =0
Liquid thermal conductivity [W/(m K)].
scalar dipm() const
Dipole moment [].
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
scalar limit(const scalar T) const
Limit the temperature to be in the range Tlow_ to Thigh_.
scalar Tc() const
Critical temperature [K].
scalar Zc() const
Critical compressibilty factor.
scalar Hc() const
Chemical enthalpy [J/kg].
scalar psi(scalar p, scalar T) const
Liquid compressibility rho/p [s^2/m^2].
A class for handling words, derived from string.
Definition: word.H:59
virtual scalar h(scalar p, scalar T) const =0
Liquid enthalpy [J/kg] - reference to 298.15 K.
virtual scalar Cp(const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/(kg K)].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
scalar Vc() const
Critical volume [m^3/kmol].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
scalar Hs(const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const volScalarField & T
virtual void writeData(Ostream &os) const =0
Write the function coefficients.
scalar omega() const
Pitzer&#39;s ascentric factor [].
scalar Pc() const
Critical pressure [Pa].
volScalarField & p
scalar alphah(const scalar p, const scalar T) const
Liquid thermal diffusivity of enthalpy [kg/ms].
scalar Tb() const
Normal boiling temperature [K].
scalar Tt() const
Triple point temperature [K].