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