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-2018 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>::S(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>::psi(scalar p, scalar T) const
126 {
127  return psi_;
128 }
129 
130 
131 template<class Specie>
132 inline Foam::scalar Foam::linear<Specie>::Z(scalar p, scalar T) const
133 {
134  return 1;
135 }
136 
137 
138 template<class Specie>
139 inline Foam::scalar Foam::linear<Specie>::CpMCv(scalar p, scalar T) const
140 {
141  return 0;
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
146 
147 template<class Specie>
148 inline void Foam::linear<Specie>::operator+=
149 (
150  const linear<Specie>& pf
151 )
152 {
153  scalar Y1 = this->Y();
154  Specie::operator+=(pf);
155 
156  if (mag(this->Y()) > small)
157  {
158  Y1 /= this->Y();
159  const scalar Y2 = pf.Y()/this->Y();
160 
161  psi_ = Y1*psi_ + Y2*pf.psi_;
162  rho0_ = Y1*rho0_ + Y2*pf.rho0_;
163  }
164 }
165 
166 
167 template<class Specie>
168 inline void Foam::linear<Specie>::operator*=(const scalar s)
169 {
170  Specie::operator*=(s);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
175 
176 template<class Specie>
177 inline Foam::linear<Specie> Foam::operator+
178 (
179  const linear<Specie>& pf1,
180  const linear<Specie>& pf2
181 )
182 {
183  Specie sp
184  (
185  static_cast<const Specie&>(pf1)
186  + static_cast<const Specie&>(pf2)
187  );
188 
189  if (mag(sp.Y()) < small)
190  {
191  return linear<Specie>
192  (
193  sp,
194  pf1.psi_,
195  pf1.rho0_
196  );
197  }
198  else
199  {
200  const scalar Y1 = pf1.Y()/sp.Y();
201  const scalar Y2 = pf2.Y()/sp.Y();
202 
203  return linear<Specie>
204  (
205  sp,
206  Y1*pf1.psi_ + Y2*pf2.psi_,
207  Y1*pf1.rho0_ + Y2*pf2.rho0_
208  );
209  }
210 }
211 
212 
213 template<class Specie>
214 inline Foam::linear<Specie> Foam::operator*
215 (
216  const scalar s,
217  const linear<Specie>& pf
218 )
219 {
220  return linear<Specie>
221  (
222  s*static_cast<const Specie&>(pf),
223  pf.psi_,
224  pf.rho0_
225  );
226 }
227 
228 
229 template<class Specie>
230 inline Foam::linear<Specie> Foam::operator==
231 (
232  const linear<Specie>& pf1,
233  const linear<Specie>& pf2
234 )
235 {
236  Specie sp
237  (
238  static_cast<const Specie&>(pf1)
239  == static_cast<const Specie&>(pf2)
240  );
241 
242  const scalar Y1 = pf1.Y()/sp.Y();
243  const scalar Y2 = pf2.Y()/sp.Y();
244 
245  return linear<Specie>
246  (
247  sp,
248  Y2*pf2.psi_ - Y1*pf1.psi_,
249  Y2*pf2.rho0_ - Y1*pf1.rho0_
250  );
251 }
252 
253 
254 // ************************************************************************* //
scalar H(const scalar p, const scalar T) const
Return enthalpy departure [J/kg].
Definition: linearI.H:88
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:158
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 rho/p [s^2/m^2].
Definition: linearI.H:125
scalar S(const scalar p, const scalar T) const
Return entropy [J/kg/K].
Definition: linearI.H:118
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 departure [J/(kg K].
Definition: linearI.H:111
void operator*=(const scalar)
Definition: linearI.H:168
const dimensionedScalar Pstd
Standard pressure.
linear(const fvMesh &mesh)
Construct from mesh.
Definition: linear.H:64
const volScalarField & T
scalar rho(scalar p, scalar T) const
Return density [kg/m^3].
Definition: linearI.H:81
PtrList< volScalarField > & Y
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: linearI.H:139
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
autoPtr< linear > clone() const
Construct and return a clone.
Definition: linearI.H:61
volScalarField & p
scalar Cp(scalar p, scalar T) const
Return Cp departure [J/(kg K].
Definition: linearI.H:95
scalar E(const scalar p, const scalar T) const
Return internal energy departure [J/kg].
Definition: linearI.H:102
scalar Z(scalar p, scalar T) const
Return compression factor [].
Definition: linearI.H:132