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-2023 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class Thermo>
80 inline Foam::scalar Foam::WLFTransport<Thermo>::mu
81 (
82  const scalar p,
83  const scalar T
84 ) const
85 {
86  return mu0_*exp(-C1_*(T - Tr_)/(C2_ + T - Tr_));
87 }
88 
89 
90 template<class Thermo>
92 (
93  const scalar p,
94  const scalar T
95 ) const
96 {
97  return this->Cp(p, T)*mu(p, T)*rPr_;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
102 
103 template<class Thermo>
105 (
106  const WLFTransport<Thermo>& wlft
107 )
108 {
109  scalar Y1 = this->Y();
110 
111  Thermo::operator+=(wlft);
112 
113  if (mag(this->Y()) > small)
114  {
115  Y1 /= this->Y();
116  scalar Y2 = wlft.Y()/this->Y();
117 
118  mu0_ = Y1*mu0_ + Y2*wlft.mu0_;
119  Tr_ = Y1*Tr_ + Y2*wlft.Tr_;
120  C1_ = Y1*C1_ + Y2*wlft.C1_;
121  C2_ = Y1*C2_ + Y2*wlft.C2_;
122  rPr_ = 1.0/(Y1/rPr_ + Y2/wlft.rPr_);
123  }
124 }
125 
126 
127 template<class Thermo>
129 (
130  const scalar s
131 )
132 {
133  Thermo::operator*=(s);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
138 
139 template<class Thermo>
140 inline Foam::WLFTransport<Thermo> Foam::operator+
141 (
142  const WLFTransport<Thermo>& wlft1,
143  const WLFTransport<Thermo>& wlft2
144 )
145 {
146  Thermo t
147  (
148  static_cast<const Thermo&>(wlft1) + static_cast<const Thermo&>(wlft2)
149  );
150 
151  if (mag(t.Y()) < small)
152  {
153  return WLFTransport<Thermo>
154  (
155  t,
156  0,
157  wlft1.mu0_,
158  wlft1.Tr_,
159  wlft1.C1_,
160  wlft1.C2_,
161  wlft1.rPr_
162  );
163  }
164  else
165  {
166  scalar Y1 = wlft1.Y()/t.Y();
167  scalar Y2 = wlft2.Y()/t.Y();
168 
169  return WLFTransport<Thermo>
170  (
171  t,
172  Y1*wlft1.mu0_ + Y2*wlft2.mu0_,
173  Y1*wlft1.Tr_ + Y2*wlft2.Tr_,
174  Y1*wlft1.C1_ + Y2*wlft2.C1_,
175  Y1*wlft1.C2_ + Y2*wlft2.C2_,
176  1.0/(Y1/wlft1.rPr_ + Y2/wlft2.rPr_)
177  );
178  }
179 }
180 
181 
182 template<class Thermo>
183 inline Foam::WLFTransport<Thermo> Foam::operator*
184 (
185  const scalar s,
186  const WLFTransport<Thermo>& wlft
187 )
188 {
189  return WLFTransport<Thermo>
190  (
191  s*static_cast<const Thermo&>(wlft),
192  wlft.mu0_,
193  wlft.Tr_,
194  wlft.C1_,
195  wlft.C2_,
196  1.0/wlft.rPr_
197  );
198 }
199 
200 
201 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
Transport package using the Williams-Landel-Ferry model for viscosity of polymer melts:
Definition: WLFTransport.H:132
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
Definition: WLFTransportI.H:81
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
Definition: WLFTransportI.H:92
WLFTransport(const Thermo &t, const scalar mu0, const scalar Tr, const scalar C1, const scalar C2, const scalar Pr)
Construct from components.
Definition: WLFTransportI.H:32
autoPtr< WLFTransport > clone() const
Construct and return a clone.
Definition: WLFTransportI.H:68
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A class for handling words, derived from string.
Definition: word.H:62
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p
PtrList< volScalarField > & Y