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