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 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
116 
117 template<class Thermo>
118 inline void Foam::WLFTransport<Thermo>::operator+=
119 (
120  const WLFTransport<Thermo>& wlft
121 )
122 {
123  scalar Y1 = this->Y();
124 
125  Thermo::operator+=(wlft);
126 
127  if (mag(this->Y()) > small)
128  {
129  Y1 /= this->Y();
130  scalar Y2 = wlft.Y()/this->Y();
131 
132  mu0_ = Y1*mu0_ + Y2*wlft.mu0_;
133  Tr_ = Y1*Tr_ + Y2*wlft.Tr_;
134  C1_ = Y1*C1_ + Y2*wlft.C1_;
135  C2_ = Y1*C2_ + Y2*wlft.C2_;
136  rPr_ = 1.0/(Y1/rPr_ + Y2/wlft.rPr_);
137  }
138 }
139 
140 
141 template<class Thermo>
142 inline void Foam::WLFTransport<Thermo>::operator*=
143 (
144  const scalar s
145 )
146 {
147  Thermo::operator*=(s);
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
152 
153 template<class Thermo>
154 inline Foam::WLFTransport<Thermo> Foam::operator+
155 (
156  const WLFTransport<Thermo>& wlft1,
157  const WLFTransport<Thermo>& wlft2
158 )
159 {
160  Thermo t
161  (
162  static_cast<const Thermo&>(wlft1) + static_cast<const Thermo&>(wlft2)
163  );
164 
165  if (mag(t.Y()) < small)
166  {
167  return WLFTransport<Thermo>
168  (
169  t,
170  0,
171  wlft1.mu0_,
172  wlft1.Tr_,
173  wlft1.C1_,
174  wlft1.C2_,
175  wlft1.rPr_
176  );
177  }
178  else
179  {
180  scalar Y1 = wlft1.Y()/t.Y();
181  scalar Y2 = wlft2.Y()/t.Y();
182 
183  return WLFTransport<Thermo>
184  (
185  t,
186  Y1*wlft1.mu0_ + Y2*wlft2.mu0_,
187  Y1*wlft1.Tr_ + Y2*wlft2.Tr_,
188  Y1*wlft1.C1_ + Y2*wlft2.C1_,
189  Y1*wlft1.C2_ + Y2*wlft2.C2_,
190  1.0/(Y1/wlft1.rPr_ + Y2/wlft2.rPr_)
191  );
192  }
193 }
194 
195 
196 template<class Thermo>
197 inline Foam::WLFTransport<Thermo> Foam::operator*
198 (
199  const scalar s,
200  const WLFTransport<Thermo>& wlft
201 )
202 {
203  return WLFTransport<Thermo>
204  (
205  s*static_cast<const Thermo&>(wlft),
206  wlft.mu0_,
207  wlft.Tr_,
208  wlft.C1_,
209  wlft.C2_,
210  1.0/wlft.rPr_
211  );
212 }
213 
214 
215 // ************************************************************************* //
dictionary dict
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
Definition: WLFTransportI.H:95
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
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].