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-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 "linear.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Specie>
48 (
49  const word& name,
50  const linear<Specie>& pf
51 )
52 :
53  Specie(name, pf),
54  psi_(pf.psi_),
55  rho0_(pf.rho0_)
56 {}
57 
58 
59 template<class Specie>
62 {
63  return autoPtr<linear<Specie>>(new linear<Specie>(*this));
64 }
65 
66 
67 template<class Specie>
70 (
71  const dictionary& dict
72 )
73 {
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class Specie>
81 inline Foam::scalar Foam::linear<Specie>::rho(scalar p, scalar T) const
82 {
83  return rho0_ + psi_*p;
84 }
85 
86 
87 template<class Specie>
88 inline Foam::scalar Foam::linear<Specie>::H(scalar p, scalar T) const
89 {
90  return log(rho(p, T))/psi_;
91 }
92 
93 
94 template<class Specie>
95 inline Foam::scalar Foam::linear<Specie>::Cp(scalar p, scalar T) const
96 {
97  return 0;
98 }
99 
100 
101 template<class Specie>
102 inline Foam::scalar Foam::linear<Specie>::E(scalar p, scalar T) const
103 {
104  const scalar rho = this->rho(p, T);
105 
106  return log(rho)/psi_ - p/rho;
107 }
108 
109 
110 template<class Specie>
111 inline Foam::scalar Foam::linear<Specie>::Cv(scalar p, scalar T) const
112 {
113  return 0;
114 }
115 
116 
117 template<class Specie>
118 inline Foam::scalar Foam::linear<Specie>::Sp(scalar p, scalar T) const
119 {
120  return -log((rho0_ + psi_*p)/(rho0_ + psi_*Pstd))/(T*psi_);
121 }
122 
123 
124 template<class Specie>
125 inline Foam::scalar Foam::linear<Specie>::Sv(scalar p, scalar T) const
126 {
128  return 0;
129 }
130 
131 
132 template<class Specie>
133 inline Foam::scalar Foam::linear<Specie>::psi(scalar p, scalar T) const
134 {
135  return psi_;
136 }
137 
138 
139 template<class Specie>
140 inline Foam::scalar Foam::linear<Specie>::Z(scalar p, scalar T) const
141 {
142  return p/(rho(p, T)*this->R()*T);
143 }
144 
145 
146 template<class Specie>
147 inline Foam::scalar Foam::linear<Specie>::CpMCv(scalar p, scalar T) const
148 {
149  return 0;
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
155 template<class Specie>
156 inline void Foam::linear<Specie>::operator+=
157 (
158  const linear<Specie>& pf
159 )
160 {
162 }
163 
164 
165 template<class Specie>
166 inline void Foam::linear<Specie>::operator*=(const scalar s)
167 {
168  Specie::operator*=(s);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
173 
174 template<class Specie>
175 inline Foam::linear<Specie> Foam::operator+
176 (
177  const linear<Specie>& pf1,
178  const linear<Specie>& pf2
179 )
180 {
182  return pf1;
183 }
184 
185 
186 template<class Specie>
187 inline Foam::linear<Specie> Foam::operator*
188 (
189  const scalar s,
190  const linear<Specie>& pf
191 )
192 {
193  return linear<Specie>
194  (
195  s*static_cast<const Specie&>(pf),
196  pf.psi_,
197  pf.rho0_
198  );
199 }
200 
201 
202 template<class Specie>
203 inline Foam::linear<Specie> Foam::operator==
204 (
205  const linear<Specie>& pf1,
206  const linear<Specie>& pf2
207 )
208 {
210  return pf1;
211 }
212 
213 
214 // ************************************************************************* //
scalar H(const scalar p, const scalar T) const
Return enthalpy contribution [J/kg].
Definition: linearI.H:88
scalar Sv(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cv/T [J/kg/K].
Definition: linearI.H:125
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static autoPtr< linear > New(const dictionary &dict)
Definition: linearI.H:70
Central-differencing interpolation scheme class.
Definition: linear.H:50
scalar psi(scalar p, scalar T) const
Return compressibility [s^2/m^2].
Definition: linearI.H:133
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: linearI.H:111
void operator*=(const scalar)
Definition: linearI.H:166
linear(const fvMesh &mesh)
Construct from mesh.
Definition: linear.H:64
const dimensionedScalar Pstd
Standard pressure.
#define noCoefficientMixing(Type)
Definition: specie.H:159
const volScalarField & T
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: linearI.H:81
#define R(A, B, C, D, E, F, K, M)
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: linearI.H:147
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
autoPtr< linear > clone() const
Construct and return a clone.
Definition: linearI.H:61
volScalarField & p
scalar Cp(scalar p, scalar T) const
Return Cp contribution [J/(kg K].
Definition: linearI.H:95
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
scalar E(const scalar p, const scalar T) const
Return internal energy contribution [J/kg].
Definition: linearI.H:102
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: linearI.H:140
scalar Sp(const scalar p, const scalar T) const
Return entropy contribution to the integral of Cp/T [J/kg/K].
Definition: linearI.H:118