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-2020 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 Thermo& t,
34  const scalar mu0,
35  const scalar Tr,
36  const scalar C1,
37  const scalar C2,
38  const scalar Pr
39 )
40 :
41  Thermo(t),
42  mu0_(mu0),
43  Tr_(Tr),
44  C1_(C1),
45  C2_(C2),
46  rPr_(scalar(1)/Pr)
47 {}
48 
49 
50 template<class Thermo>
52 (
53  const word& name,
54  const WLFTransport& wlft
55 )
56 :
57  Thermo(name, wlft),
58  mu0_(wlft.mu0_),
59  Tr_(wlft.Tr_),
60  C1_(wlft.C1_),
61  C2_(wlft.C2_),
62  rPr_(wlft.rPr_)
63 {}
64 
65 
66 template<class Thermo>
69 {
71  (
72  new WLFTransport<Thermo>(*this)
73  );
74 }
75 
76 
77 template<class Thermo>
80 (
81  const dictionary& dict
82 )
83 {
85  (
87  );
88 }
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class Thermo>
94 inline Foam::scalar Foam::WLFTransport<Thermo>::mu
95 (
96  const scalar p,
97  const scalar T
98 ) const
99 {
100  return mu0_*exp(-C1_*(T - Tr_)/(C2_ + T - Tr_));
101 }
102 
103 
104 template<class Thermo>
105 inline Foam::scalar Foam::WLFTransport<Thermo>::kappa
106 (
107  const scalar p,
108  const scalar T
109 ) const
110 {
111  return this->Cp(p, T)*mu(p, T)*rPr_;
112 }
113 
114 
115 template<class Thermo>
116 inline Foam::scalar Foam::WLFTransport<Thermo>::alphah
117 (
118  const scalar p,
119  const scalar T
120 ) const
121 {
122 
123  return mu(p, T)*rPr_;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
128 
129 template<class Thermo>
130 inline void Foam::WLFTransport<Thermo>::operator+=
131 (
132  const WLFTransport<Thermo>& wlft
133 )
134 {
135  scalar Y1 = this->Y();
136 
137  Thermo::operator+=(wlft);
138 
139  if (mag(this->Y()) > small)
140  {
141  Y1 /= this->Y();
142  scalar Y2 = wlft.Y()/this->Y();
143 
144  mu0_ = Y1*mu0_ + Y2*wlft.mu0_;
145  Tr_ = Y1*Tr_ + Y2*wlft.Tr_;
146  C1_ = Y1*C1_ + Y2*wlft.C1_;
147  C2_ = Y1*C2_ + Y2*wlft.C2_;
148  rPr_ = 1.0/(Y1/rPr_ + Y2/wlft.rPr_);
149  }
150 }
151 
152 
153 template<class Thermo>
154 inline void Foam::WLFTransport<Thermo>::operator*=
155 (
156  const scalar s
157 )
158 {
159  Thermo::operator*=(s);
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
164 
165 template<class Thermo>
166 inline Foam::WLFTransport<Thermo> Foam::operator+
167 (
168  const WLFTransport<Thermo>& wlft1,
169  const WLFTransport<Thermo>& wlft2
170 )
171 {
172  Thermo t
173  (
174  static_cast<const Thermo&>(wlft1) + static_cast<const Thermo&>(wlft2)
175  );
176 
177  if (mag(t.Y()) < small)
178  {
179  return WLFTransport<Thermo>
180  (
181  t,
182  0,
183  wlft1.mu0_,
184  wlft1.Tr_,
185  wlft1.C1_,
186  wlft1.C2_,
187  wlft1.rPr_
188  );
189  }
190  else
191  {
192  scalar Y1 = wlft1.Y()/t.Y();
193  scalar Y2 = wlft2.Y()/t.Y();
194 
195  return WLFTransport<Thermo>
196  (
197  t,
198  Y1*wlft1.mu0_ + Y2*wlft2.mu0_,
199  Y1*wlft1.Tr_ + Y2*wlft2.Tr_,
200  Y1*wlft1.C1_ + Y2*wlft2.C1_,
201  Y1*wlft1.C2_ + Y2*wlft2.C2_,
202  1.0/(Y1/wlft1.rPr_ + Y2/wlft2.rPr_)
203  );
204  }
205 }
206 
207 
208 template<class Thermo>
209 inline Foam::WLFTransport<Thermo> Foam::operator*
210 (
211  const scalar s,
212  const WLFTransport<Thermo>& wlft
213 )
214 {
215  return WLFTransport<Thermo>
216  (
217  s*static_cast<const Thermo&>(wlft),
218  wlft.mu0_,
219  wlft.Tr_,
220  wlft.C1_,
221  wlft.C2_,
222  1.0/wlft.rPr_
223  );
224 }
225 
226 
227 // ************************************************************************* //
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:95
const dimensionedScalar & mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
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].
autoPtr< WLFTransport > clone() const
Construct and return a clone.
Definition: WLFTransportI.H:68
static autoPtr< WLFTransport > New(const dictionary &dict)
Definition: WLFTransportI.H:80
Transport package using the Williams-Landel-Ferry model.
Definition: WLFTransport.H:61
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].