icoTabulatedI.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) 2020-2023 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 #include "icoTabulated.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const nonUniformTable& rho
35 )
36 :
37  Specie(sp),
38  rho_(rho)
39 {}
40 
41 
42 template<class Specie>
44 (
45  const word& name,
46  const icoTabulated<Specie>& ip
47 )
48 :
49  Specie(name, ip),
50  rho_(ip.rho_)
51 {}
52 
53 
54 template<class Specie>
57 {
59  (
60  new icoTabulated<Specie>(*this)
61  );
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class Specie>
69 (
70  scalar p,
71  scalar T
72 ) const
73 {
74  return rho_.value(T);
75 }
76 
77 
78 template<class Specie>
79 inline Foam::scalar Foam::icoTabulated<Specie>::H
80 (
81  scalar p,
82  scalar T
83 ) const
84 {
85  return p/this->rho(p, T);
86 }
87 
88 
89 template<class Specie>
90 inline Foam::scalar Foam::icoTabulated<Specie>::Cp
91 (
92  scalar p,
93  scalar T
94 ) const
95 {
96  return 0;
97 }
98 
99 
100 template<class Specie>
101 inline Foam::scalar Foam::icoTabulated<Specie>::E
102 (
103  scalar p,
104  scalar T
105 ) const
106 {
107  return 0;
108 }
109 
110 
111 template<class Specie>
113 (
114  scalar p,
115  scalar T
116 ) const
117 {
118  return 0;
119 }
120 
121 
122 template<class Specie>
124 (
125  scalar p,
126  scalar T
127 ) const
128 {
129  return 0;
130 }
131 
132 
133 template<class Specie>
135 (
136  scalar p,
137  scalar T
138 ) const
139 {
140  return 0;
141 }
142 
143 
144 template<class Specie>
146 (
147  scalar p,
148  scalar T
149 ) const
150 {
151  return 0;
152 }
153 
154 
155 template<class Specie>
156 inline Foam::scalar Foam::icoTabulated<Specie>::Z
157 (
158  scalar p,
159  scalar T
160 ) const
161 {
162  return p/(rho(p, T)*this->R()*T);
163 }
164 
165 
166 template<class Specie>
168 (
169  scalar p,
170  scalar T
171 ) const
172 {
173  return 0;
174 }
175 
176 
177 template<class Specie>
179 (
180  scalar p,
181  scalar T
182 ) const
183 {
184  return -rho_.dfdT(T)/rho(p, T);
185 }
186 
187 
188 // ************************************************************************* //
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Incompressible of equation of state using non-uniform tabulated density vs temperature.
Definition: icoTabulated.H:98
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: icoTabulatedI.H:80
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
icoTabulated(const Specie &sp, const nonUniformTable &rho)
Construct from components.
Definition: icoTabulatedI.H:32
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: icoTabulatedI.H:69
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
autoPtr< icoTabulated > clone() const
Construct and return a clone.
Definition: icoTabulatedI.H:56
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: icoTabulatedI.H:91
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
scalar Z(scalar p, scalar T) const
Return compression factor [].
A class for handling words, derived from string.
Definition: word.H:62
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p