rhoTabulatedI.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-2021 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 "rhoTabulated.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const table2D& rho
35 )
36 :
37  Specie(sp),
38  rho_(rho)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class Specie>
46 (
47  const word& name,
48  const rhoTabulated<Specie>& ip
49 )
50 :
51  Specie(name, ip),
52  rho_(ip.rho_)
53 {}
54 
55 
56 template<class Specie>
59 {
61  (
62  new rhoTabulated<Specie>(*this)
63  );
64 }
65 
66 
67 template<class Specie>
70 {
72  (
74  );
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class Specie>
81 inline Foam::scalar Foam::rhoTabulated<Specie>::rho
82 (
83  scalar p,
84  scalar T
85 ) const
86 {
87  return rho_.value(p, T);
88 }
89 
90 
91 template<class Specie>
92 inline Foam::scalar Foam::rhoTabulated<Specie>::H
93 (
94  scalar p,
95  scalar T
96 ) const
97 {
99  return 0;
100 }
101 
102 
103 template<class Specie>
104 inline Foam::scalar Foam::rhoTabulated<Specie>::Cp
105 (
106  scalar p,
107  scalar T
108 ) const
109 {
111  return 0;
112 }
113 
114 
115 template<class Specie>
116 inline Foam::scalar Foam::rhoTabulated<Specie>::E
117 (
118  scalar p,
119  scalar T
120 ) const
121 {
123  return 0;
124 }
125 
126 
127 template<class Specie>
128 inline Foam::scalar Foam::rhoTabulated<Specie>::Cv
129 (
130  scalar p,
131  scalar T
132 ) const
133 {
135  return 0;
136 }
137 
138 
139 template<class Specie>
140 inline Foam::scalar Foam::rhoTabulated<Specie>::Sp
141 (
142  scalar p,
143  scalar T
144 ) const
145 {
147  return 0;
148 }
149 
150 
151 template<class Specie>
152 inline Foam::scalar Foam::rhoTabulated<Specie>::Sv
153 (
154  scalar p,
155  scalar T
156 ) const
157 {
159  return 0;
160 }
161 
162 
163 template<class Specie>
164 inline Foam::scalar Foam::rhoTabulated<Specie>::psi
165 (
166  scalar p,
167  scalar T
168 ) const
169 {
170  return rho_.dfdp(p, T);
171 }
172 
173 
174 template<class Specie>
175 inline Foam::scalar Foam::rhoTabulated<Specie>::Z
176 (
177  scalar p,
178  scalar T
179 ) const
180 {
181  return p/(rho(p, T)*this->R()*T);
182 }
183 
184 
185 template<class Specie>
186 inline Foam::scalar Foam::rhoTabulated<Specie>::CpMCv
187 (
188  scalar p,
189  scalar T
190 ) const
191 {
193  return 0;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
198 
199 template<class Specie>
200 inline void Foam::rhoTabulated<Specie>::operator+=
201 (
202  const rhoTabulated<Specie>& pf
203 )
204 {
206 }
207 
208 
209 template<class Specie>
210 inline void Foam::rhoTabulated<Specie>::operator*=(const scalar s)
211 {
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
217 
218 template<class Specie>
219 inline Foam::rhoTabulated<Specie> Foam::operator+
220 (
221  const rhoTabulated<Specie>& pf1,
222  const rhoTabulated<Specie>& pf2
223 )
224 {
226  return pf1;
227 }
228 
229 
230 template<class Specie>
231 inline Foam::rhoTabulated<Specie> Foam::operator*
232 (
233  const scalar s,
234  const rhoTabulated<Specie>& pf
235 )
236 {
238  return pf;
239 }
240 
241 
242 template<class Specie>
243 inline Foam::rhoTabulated<Specie> Foam::operator==
244 (
245  const rhoTabulated<Specie>& pf1,
246  const rhoTabulated<Specie>& pf2
247 )
248 {
250  return pf1;
251 }
252 
253 
254 // ************************************************************************* //
dictionary dict
scalar Z(scalar p, scalar T) const
Return compression factor [].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: rhoTabulatedI.H:82
static autoPtr< rhoTabulated > New(const dictionary &dict)
Definition: rhoTabulatedI.H:69
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
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: rhoTabulatedI.H:93
A class for handling words, derived from string.
Definition: word.H:59
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Incompressible of equation of state using uniform tabulated density vs pressure and temperature...
Definition: rhoTabulated.H:96
autoPtr< rhoTabulated > clone() const
Construct and return a clone.
Definition: rhoTabulatedI.H:58
#define noCoefficientMixing(Type)
Definition: specie.H:159
const volScalarField & T
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
#define R(A, B, C, D, E, F, K, M)
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
rhoTabulated(const Specie &sp, const table2D &rho)
Construct from components.
Definition: rhoTabulatedI.H:32
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
void operator*=(const scalar)