WLFTransportI.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) 2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "specie.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Thermo>
34 (
35  const word& name,
36  const WLFTransport& wlft
37 )
38 :
39  Thermo(name, wlft),
40  mu0_(wlft.mu0_),
41  Tr_(wlft.Tr_),
42  C1_(wlft.C1_),
43  C2_(wlft.C2_),
44  rPr_(wlft.rPr_)
45 {}
46 
47 
48 template<class Thermo>
51 {
53  (
54  new WLFTransport<Thermo>(*this)
55  );
56 }
57 
58 
59 template<class Thermo>
62 (
63  const dictionary& dict
64 )
65 {
67  (
69  );
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class Thermo>
76 inline Foam::scalar Foam::WLFTransport<Thermo>::mu
77 (
78  const scalar p,
79  const scalar T
80 ) const
81 {
82  return mu0_*exp(-C1_*(T - Tr_)/(C2_ + T - Tr_));
83 }
84 
85 
86 template<class Thermo>
87 inline Foam::scalar Foam::WLFTransport<Thermo>::kappa
88 (
89  const scalar p, const scalar T
90 ) const
91 {
92  return this->Cp(p, T)*mu(p, T)*rPr_;
93 }
94 
95 
96 template<class Thermo>
97 inline Foam::scalar Foam::WLFTransport<Thermo>::alphah
98 (
99  const scalar p,
100  const scalar T
101 ) const
102 {
103 
104  return mu(p, T)*rPr_;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
109 
110 template<class Thermo>
112 (
113  const WLFTransport<Thermo>& wlft
114 )
115 {
116  scalar Y1 = this->Y();
117 
118  Thermo::operator+=(wlft);
119 
120  if (mag(this->Y()) > SMALL)
121  {
122  Y1 /= this->Y();
123  scalar Y2 = wlft.Y()/this->Y();
124 
125  mu0_ = Y1*mu0_ + Y2*wlft.mu0_;
126  Tr_ = Y1*Tr_ + Y2*wlft.Tr_;
127  C1_ = Y1*C1_ + Y2*wlft.C1_;
128  C2_ = Y1*C2_ + Y2*wlft.C2_;
129  rPr_ = 1.0/(Y1/rPr_ + Y2/wlft.rPr_);
130  }
131 }
132 
133 
134 template<class Thermo>
136 (
137  const scalar s
138 )
139 {
140  Thermo::operator*=(s);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
145 
146 template<class Thermo>
147 inline Foam::WLFTransport<Thermo> Foam::operator+
148 (
149  const WLFTransport<Thermo>& wlft1,
150  const WLFTransport<Thermo>& wlft2
151 )
152 {
153  Thermo t
154  (
155  static_cast<const Thermo&>(wlft1) + static_cast<const Thermo&>(wlft2)
156  );
157 
158  if (mag(t.Y()) < SMALL)
159  {
160  return WLFTransport<Thermo>
161  (
162  t,
163  0,
164  wlft1.mu0_,
165  wlft1.Tr_,
166  wlft1.C1_,
167  wlft1.C2_,
168  wlft1.rPr_
169  );
170  }
171  else
172  {
173  scalar Y1 = wlft1.Y()/t.Y();
174  scalar Y2 = wlft2.Y()/t.Y();
175 
176  return WLFTransport<Thermo>
177  (
178  t,
179  Y1*wlft1.mu0_ + Y2*wlft2.mu0_,
180  Y1*wlft1.Tr_ + Y2*wlft2.Tr_,
181  Y1*wlft1.C1_ + Y2*wlft2.C1_,
182  Y1*wlft1.C2_ + Y2*wlft2.C2_,
183  1.0/(Y1/wlft1.rPr_ + Y2/wlft2.rPr_)
184  );
185  }
186 }
187 
188 
189 template<class Thermo>
190 inline Foam::WLFTransport<Thermo> Foam::operator*
191 (
192  const scalar s,
193  const WLFTransport<Thermo>& wlft
194 )
195 {
196  return WLFTransport<Thermo>
197  (
198  s*static_cast<const Thermo&>(wlft),
199  wlft.mu0_,
200  wlft.Tr_,
201  wlft.C1_,
202  wlft.C2_,
203  1.0/wlft.rPr_
204  );
205 }
206 
207 
208 // ************************************************************************* //
Foam::WLFTransport::alphah
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
Definition: WLFTransportI.H:98
Foam::WLFTransport::New
static autoPtr< WLFTransport > New(const dictionary &dict)
Definition: WLFTransportI.H:62
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
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))
Definition: gmvOutputSpray.H:25
Foam::WLFTransport::clone
autoPtr< WLFTransport > clone() const
Construct and return a clone.
Definition: WLFTransportI.H:50
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
specie.H
Foam::WLFTransport::WLFTransport
WLFTransport(const word &, const WLFTransport &)
Construct as named copy.
Definition: WLFTransportI.H:34
Foam::WLFTransport
Transport package using the Williams-Landel-Ferry model.
Definition: WLFTransport.H:64
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::WLFTransport::mu
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
Definition: WLFTransportI.H:77
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::WLFTransport::kappa
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
Definition: WLFTransportI.H:88