linearI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "linear.H"
29
30// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31
32template<class Specie>
34(
35 const Specie& sp,
36 const scalar psi,
37 const scalar rho0
38)
39:
40 Specie(sp),
41 psi_(psi),
42 rho0_(rho0)
43{}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48template<class Specie>
50(
51 const word& name,
52 const linear<Specie>& pf
53)
54:
55 Specie(name, pf),
56 psi_(pf.psi_),
57 rho0_(pf.rho0_)
58{}
59
60
61template<class Specie>
64{
65 return autoPtr<linear<Specie>>::New(*this);
66}
67
68
69template<class Specie>
72(
73 const dictionary& dict
74)
75{
77}
78
79
80// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81
82template<class Specie>
83inline Foam::scalar Foam::linear<Specie>::rho(scalar p, scalar T) const
84{
85 return rho0_ + psi_*p;
86}
87
88
89template<class Specie>
90inline Foam::scalar Foam::linear<Specie>::H(scalar p, scalar T) const
91{
92 return 0;
93}
94
95
96template<class Specie>
97inline Foam::scalar Foam::linear<Specie>::Cp(scalar p, scalar T) const
98{
99 return 0;
100}
101
102template<class Specie>
103inline Foam::scalar Foam::linear<Specie>::E(scalar p, scalar T) const
104{
105 return 0;
106}
107
108
109template<class Specie>
110inline Foam::scalar Foam::linear<Specie>::Cv(scalar p, scalar T) const
111{
112 return 0;
113}
114
115
116template<class Specie>
117inline Foam::scalar Foam::linear<Specie>::S(scalar p, scalar T) const
118{
119 return -log((rho0_ + psi_*p)/(rho0_ + psi_*Pstd))/(T*psi_);
120}
121
122
123template<class Specie>
124inline Foam::scalar Foam::linear<Specie>::psi(scalar p, scalar T) const
125{
126 return psi_;
127}
128
129
130template<class Specie>
131inline Foam::scalar Foam::linear<Specie>::Z(scalar p, scalar T) const
132{
133 return 1;
134}
135
136
137template<class Specie>
138inline Foam::scalar Foam::linear<Specie>::CpMCv(scalar p, scalar T) const
139{
140 return 0;
141}
142
143
144// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
145
146template<class Specie>
148(
149 const linear<Specie>& pf
150)
151{
152 scalar Y1 = this->Y();
153 Specie::operator+=(pf);
154
155 if (mag(this->Y()) > SMALL)
156 {
157 Y1 /= this->Y();
158 const scalar Y2 = pf.Y()/this->Y();
159
160 psi_ = Y1*psi_ + Y2*pf.psi_;
161 rho0_ = Y1*rho0_ + Y2*pf.rho0_;
162 }
163}
164
165
166template<class Specie>
167inline void Foam::linear<Specie>::operator*=(const scalar s)
168{
169 Specie::operator*=(s);
170}
171
172
173// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
174
175template<class Specie>
176inline Foam::linear<Specie> Foam::operator+
177(
178 const linear<Specie>& pf1,
179 const linear<Specie>& pf2
180)
181{
182 Specie sp
183 (
184 static_cast<const Specie&>(pf1)
185 + static_cast<const Specie&>(pf2)
186 );
187
188 if (mag(sp.Y()) < SMALL)
189 {
190 return linear<Specie>
191 (
192 sp,
193 pf1.psi_,
194 pf1.rho0_
195 );
196 }
197 else
198 {
199 const scalar Y1 = pf1.Y()/sp.Y();
200 const scalar Y2 = pf2.Y()/sp.Y();
201
202 return linear<Specie>
203 (
204 sp,
205 Y1*pf1.psi_ + Y2*pf2.psi_,
206 Y1*pf1.rho0_ + Y2*pf2.rho0_
207 );
208 }
209}
210
211
212template<class Specie>
213inline Foam::linear<Specie> Foam::operator*
214(
215 const scalar s,
216 const linear<Specie>& pf
217)
218{
219 return linear<Specie>
220 (
221 s*static_cast<const Specie&>(pf),
222 pf.psi_,
223 pf.rho0_
224 );
225}
226
227
228template<class Specie>
229inline Foam::linear<Specie> Foam::operator==
230(
231 const linear<Specie>& pf1,
232 const linear<Specie>& pf2
233)
234{
235 Specie sp
236 (
237 static_cast<const Specie&>(pf1)
238 == static_cast<const Specie&>(pf2)
239 );
240
241 const scalar Y1 = pf1.Y()/sp.Y();
242 const scalar Y2 = pf2.Y()/sp.Y();
243
244 return linear<Specie>
245 (
246 sp,
247 Y2*pf2.psi_ - Y1*pf1.psi_,
248 Y2*pf2.rho0_ - Y1*pf1.rho0_
249 );
250}
251
252
253// ************************************************************************* //
compactSpatialTensor S
The joint motion sub-space (3-DoF)
Definition: joint.H:132
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
Central-differencing interpolation scheme class.
Definition: linear.H:58
scalar CpMCv(scalar p, scalar T) const
Return (Cp - Cv) [J/(kg K].
Definition: linearI.H:138
autoPtr< linear > clone() const
Construct and return a clone.
Definition: linearI.H:63
void operator*=(const scalar)
Definition: linearI.H:167
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & psi
const volScalarField & T
const volScalarField & Cv
Definition: EEqn.H:8
const volScalarField & Cp
Definition: EEqn.H:7
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
scalar rho0
dictionary dict