constTransportI.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) 2011-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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29
30template<class Thermo>
32(
33 const Thermo& t,
34 const scalar mu,
35 const scalar Pr
36)
37:
38 Thermo(t),
39 mu_(mu),
40 rPr_(1.0/Pr)
41{}
42
43
44template<class Thermo>
46(
47 const word& name,
48 const constTransport& ct
49)
50:
51 Thermo(name, ct),
52 mu_(ct.mu_),
53 rPr_(ct.rPr_)
54{}
55
56
57template<class Thermo>
60{
62}
63
64
65template<class Thermo>
68(
69 const dictionary& dict
70)
71{
73}
74
75
76// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77
78template<class Thermo>
80(
81 const scalar p,
82 const scalar T
83) const
84{
85 return mu_;
86}
87
88
89template<class Thermo>
91(
92 const scalar p,
93 const scalar T
94) const
95{
96 return this->Cp(p, T)*mu(p, T)*rPr_;
97}
98
99
100template<class Thermo>
102(
103 const scalar p,
104 const scalar T
105) const
106{
107 return mu(p, T)*rPr_;
108}
109
110
111// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
112
113template<class Thermo>
115(
116 const constTransport<Thermo>& st
117)
118{
119 scalar Y1 = this->Y();
120
121 Thermo::operator+=(st);
122
123 if (mag(this->Y()) > SMALL)
124 {
125 Y1 /= this->Y();
126 scalar Y2 = st.Y()/this->Y();
127
128 mu_ = Y1*mu_ + Y2*st.mu_;
129 rPr_ = 1.0/(Y1/rPr_ + Y2/st.rPr_);
130 }
131}
132
133
134template<class Thermo>
136(
137 const scalar s
138)
139{
140 Thermo::operator*=(s);
141}
142
143
144// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
145
146template<class Thermo>
147inline Foam::constTransport<Thermo> Foam::operator+
148(
149 const constTransport<Thermo>& ct1,
150 const constTransport<Thermo>& ct2
151)
152{
153 Thermo t
154 (
155 static_cast<const Thermo&>(ct1) + static_cast<const Thermo&>(ct2)
156 );
157
158 if (mag(t.Y()) < SMALL)
159 {
161 (
162 t,
163 0,
164 ct1.rPr_
165 );
166 }
167 else
168 {
169 scalar Y1 = ct1.Y()/t.Y();
170 scalar Y2 = ct2.Y()/t.Y();
171
172 return constTransport<Thermo>
173 (
174 t,
175 Y1*ct1.mu_ + Y2*ct2.mu_,
176 1.0/(Y1/ct1.rPr_ + Y2/ct2.rPr_)
177 );
178 }
179}
180
181
182template<class Thermo>
183inline Foam::constTransport<Thermo> Foam::operator*
184(
185 const scalar s,
186 const constTransport<Thermo>& ct
187)
188{
189 return constTransport<Thermo>
190 (
191 s*static_cast<const Thermo&>(ct),
192 ct.mu_,
193 1.0/ct.rPr_
194 );
195}
196
197
198// ************************************************************************* //
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
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
autoPtr< constTransport > clone() const
Construct and return a clone.
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
const volScalarField & mu
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))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
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
dictionary dict
dimensionedScalar Pr("Pr", dimless, laminarTransport)