linearI.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) 2013-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 "linear.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Specie>
32 (
33  const Specie& sp,
34  const scalar psi,
35  const scalar rho0
36 )
37 :
38  Specie(sp),
39  psi_(psi),
40  rho0_(rho0)
41 {}
42 
43 
44 template<class Specie>
46 (
47  const word& name,
48  const linear<Specie>& pf
49 )
50 :
51  Specie(name, pf),
52  psi_(pf.psi_),
53  rho0_(pf.rho0_)
54 {}
55 
56 
57 template<class Specie>
60 {
61  return autoPtr<linear<Specie>>(new linear<Specie>(*this));
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class Specie>
68 inline Foam::scalar Foam::linear<Specie>::rho(scalar p, scalar T) const
69 {
70  return rho0_ + psi_*p;
71 }
72 
73 
74 template<class Specie>
75 inline Foam::scalar Foam::linear<Specie>::h(scalar p, scalar T) const
76 {
77  return log(rho(p, T))/psi_;
78 }
79 
80 
81 template<class Specie>
82 inline Foam::scalar Foam::linear<Specie>::Cp(scalar p, scalar T) const
83 {
84  return 0;
85 }
86 
87 
88 template<class Specie>
89 inline Foam::scalar Foam::linear<Specie>::e(scalar p, scalar T) const
90 {
91  const scalar rho = this->rho(p, T);
92 
93  return log(rho)/psi_ - p/rho;
94 }
95 
96 
97 template<class Specie>
98 inline Foam::scalar Foam::linear<Specie>::Cv(scalar p, scalar T) const
99 {
100  return 0;
101 }
102 
103 
104 template<class Specie>
105 inline Foam::scalar Foam::linear<Specie>::sp(scalar p, scalar T) const
106 {
107  return -log((rho0_ + psi_*p)/(rho0_ + psi_*Pstd))/(T*psi_);
108 }
109 
110 
111 template<class Specie>
112 inline Foam::scalar Foam::linear<Specie>::sv(scalar p, scalar T) const
113 {
115  return 0;
116 }
117 
118 
119 template<class Specie>
120 inline Foam::scalar Foam::linear<Specie>::psi(scalar p, scalar T) const
121 {
122  return psi_;
123 }
124 
125 
126 template<class Specie>
127 inline Foam::scalar Foam::linear<Specie>::Z(scalar p, scalar T) const
128 {
129  return p/(rho(p, T)*this->R()*T);
130 }
131 
132 
133 template<class Specie>
134 inline Foam::scalar Foam::linear<Specie>::CpMCv(scalar p, scalar T) const
135 {
136  return 0;
137 }
138 
139 
140 template<class Specie>
141 inline Foam::scalar Foam::linear<Specie>::alphav(scalar p, scalar T) const
142 {
143  return 0;
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
148 
149 template<class Specie>
151 (
152  const linear<Specie>& pf
153 )
154 {
156 }
157 
158 
159 template<class Specie>
160 inline void Foam::linear<Specie>::operator*=(const scalar s)
161 {
162  Specie::operator*=(s);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
167 
168 template<class Specie>
169 inline Foam::linear<Specie> Foam::operator+
170 (
171  const linear<Specie>& pf1,
172  const linear<Specie>& pf2
173 )
174 {
176  return pf1;
177 }
178 
179 
180 template<class Specie>
181 inline Foam::linear<Specie> Foam::operator*
182 (
183  const scalar s,
184  const linear<Specie>& pf
185 )
186 {
187  return linear<Specie>
188  (
189  s*static_cast<const Specie&>(pf),
190  pf.psi_,
191  pf.rho0_
192  );
193 }
194 
195 
196 template<class Specie>
197 inline Foam::linear<Specie> Foam::operator==
198 (
199  const linear<Specie>& pf1,
200  const linear<Specie>& pf2
201 )
202 {
203  noCoefficientMixing(linear);
204  return pf1;
205 }
206 
207 
208 // ************************************************************************* //
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Centred interpolation interpolation scheme class.
Definition: linear.H:53
scalar Cv(scalar p, scalar T) const
Return Cv contribution [J/(kg K].
Definition: linearI.H:98
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: linearI.H:120
scalar alphav(const scalar p, const scalar T) const
Return volumetric coefficient of thermal expansion [1/T].
Definition: linearI.H:141
scalar e(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: linearI.H:89
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: linearI.H:68
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: linearI.H:134
scalar h(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: linearI.H:75
linear(const fvMesh &mesh)
Construct from mesh.
Definition: linear.H:64
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: linearI.H:82
scalar sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: linearI.H:112
scalar sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: linearI.H:105
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: linearI.H:127
autoPtr< linear > clone() const
Construct and return a clone.
Definition: linearI.H:59
void operator*=(const scalar)
Definition: linearI.H:160
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
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))
const volScalarField & psi
const dimensionedScalar Pstd
Standard pressure.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedScalar log(const dimensionedScalar &ds)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar rho0
volScalarField & p
#define noCoefficientMixing(Type)
Definition: specie.H:159