rhoConstI.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) 2012-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 "rhoConst.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const scalar 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 rhoConst<Specie>& rc
49 )
50 :
51  Specie(name, rc),
52  rho_(rc.rho_)
53 {}
54 
55 
56 template<class Specie>
59 {
60  return autoPtr<rhoConst<Specie>>(new rhoConst<Specie>(*this));
61 }
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
66 template<class Specie>
67 inline Foam::scalar Foam::rhoConst<Specie>::rho(scalar p, scalar T) const
68 {
69  return rho_;
70 }
71 
72 
73 template<class Specie>
74 inline Foam::scalar Foam::rhoConst<Specie>::H(scalar p, scalar T) const
75 {
76  return p/rho_;
77 }
78 
79 
80 template<class Specie>
81 inline Foam::scalar Foam::rhoConst<Specie>::Cp(scalar p, scalar T) const
82 {
83  return 0;
84 }
85 
86 
87 template<class Specie>
88 inline Foam::scalar Foam::rhoConst<Specie>::E(scalar p, scalar T) const
89 {
90  return 0;
91 }
92 
93 
94 template<class Specie>
95 inline Foam::scalar Foam::rhoConst<Specie>::Cv(scalar p, scalar T) const
96 {
97  return 0;
98 }
99 
100 
101 template<class Specie>
102 inline Foam::scalar Foam::rhoConst<Specie>::Sp(scalar p, scalar T) const
103 {
104  return 0;
105 }
106 
107 
108 template<class Specie>
109 inline Foam::scalar Foam::rhoConst<Specie>::Sv(scalar p, scalar T) const
110 {
111  return 0;
112 }
113 
114 
115 template<class Specie>
116 inline Foam::scalar Foam::rhoConst<Specie>::psi(scalar p, scalar T) const
117 {
118  return 0;
119 }
120 
121 
122 template<class Specie>
123 inline Foam::scalar Foam::rhoConst<Specie>::Z(scalar p, scalar T) const
124 {
125  return p/(rho(p, T)*this->R()*T);
126 }
127 
128 
129 template<class Specie>
130 inline Foam::scalar Foam::rhoConst<Specie>::CpMCv(scalar p, scalar T) const
131 {
132  return 0;
133 }
134 
135 
136 template<class Specie>
137 inline Foam::scalar Foam::rhoConst<Specie>::alphav(scalar p, scalar T) const
138 {
139  return 0;
140 }
141 
142 
143 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
144 
145 template<class Specie>
147 {
148  const scalar Y1 = this->Y();
149  Specie::operator+=(rc);
150 
151  if (mag(this->Y()) > small)
152  {
153  rho_ = this->Y()/(Y1/rho_ + rc.Y()/rc.rho_);
154  }
155 }
156 
157 
158 template<class Specie>
159 inline void Foam::rhoConst<Specie>::operator*=(const scalar s)
160 {
161  Specie::operator*=(s);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
166 
167 template<class Specie>
168 inline Foam::rhoConst<Specie> Foam::operator+
169 (
170  const rhoConst<Specie>& rc1,
171  const rhoConst<Specie>& rc2
172 )
173 {
174  Specie sp
175  (
176  static_cast<const Specie&>(rc1)
177  + static_cast<const Specie&>(rc2)
178  );
179 
180  if (mag(sp.Y()) < small)
181  {
182  return rhoConst<Specie>
183  (
184  sp,
185  rc1.rho_
186  );
187  }
188  else
189  {
190  return rhoConst<Specie>
191  (
192  sp,
193  sp.Y()/(rc1.Y()/rc1.rho_ + rc2.Y()/rc2.rho_)
194  );
195  }
196 }
197 
198 
199 template<class Specie>
200 inline Foam::rhoConst<Specie> Foam::operator*
201 (
202  const scalar s,
203  const rhoConst<Specie>& rc
204 )
205 {
206  return rhoConst<Specie>(s*static_cast<const Specie&>(rc), rc.rho_);
207 }
208 
209 
210 template<class Specie>
211 inline Foam::rhoConst<Specie> Foam::operator==
212 (
213  const rhoConst<Specie>& rc1,
214  const rhoConst<Specie>& rc2
215 )
216 {
217  return rhoConst<Specie>
218  (
219  static_cast<const Specie&>(rc1) == static_cast<const Specie&>(rc2),
220  NaN
221  );
222 }
223 
224 
225 // ************************************************************************* //
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: rhoConstI.H:137
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: rhoConstI.H:102
void operator+=(const rhoConst &)
Definition: rhoConstI.H:146
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: rhoConstI.H:88
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: rhoConstI.H:116
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: rhoConstI.H:74
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: rhoConstI.H:123
rhoConst(const Specie &sp, const scalar rho)
Construct from components.
Definition: rhoConstI.H:32
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: rhoConstI.H:130
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: rhoConstI.H:81
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: rhoConstI.H:67
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))
A class for handling words, derived from string.
Definition: word.H:59
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: rhoConstI.H:95
autoPtr< rhoConst > clone() const
Construct and return a clone.
Definition: rhoConstI.H:58
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: rhoConstI.H:109
const volScalarField & T
#define R(A, B, C, D, E, F, K, M)
PtrList< volScalarField > & Y
void operator*=(const scalar)
Definition: rhoConstI.H:159
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Constant density equations of state.
Definition: rhoConst.H:68
volScalarField & p