sutherlandTransportI.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) 2011-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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Thermo>
32 (
33  const scalar mu1, const scalar T1,
34  const scalar mu2, const scalar T2
35 )
36 {
37  scalar rootT1 = sqrt(T1);
38  scalar mu1rootT2 = mu1*sqrt(T2);
39  scalar mu2rootT1 = mu2*rootT1;
40 
41  Ts_ = (mu2rootT1 - mu1rootT2)/(mu1rootT2/T1 - mu2rootT1/T2);
42 
43  As_ = mu1*(1.0 + Ts_/T1)/rootT1;
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class Thermo>
51 (
52  const Thermo& t,
53  const scalar As,
54  const scalar Ts
55 )
56 :
57  Thermo(t),
58  As_(As),
59  Ts_(Ts)
60 {}
61 
62 
63 template<class Thermo>
65 (
66  const Thermo& t,
67  const scalar mu1, const scalar T1,
68  const scalar mu2, const scalar T2
69 )
70 :
71  Thermo(t)
72 {
73  calcCoeffs(mu1, T1, mu2, T2);
74 }
75 
76 
77 template<class Thermo>
79 (
80  const word& name,
81  const sutherlandTransport& st
82 )
83 :
84  Thermo(name, st),
85  As_(st.As_),
86  Ts_(st.Ts_)
87 {}
88 
89 
90 template<class Thermo>
93 {
95  (
97  );
98 }
99 
100 
101 template<class Thermo>
104 (
105  const dictionary& dict
106 )
107 {
109  (
111  );
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class Thermo>
118 inline Foam::scalar Foam::sutherlandTransport<Thermo>::mu
119 (
120  const scalar p,
121  const scalar T
122 ) const
123 {
124  return As_*::sqrt(T)/(1.0 + Ts_/T);
125 }
126 
127 
128 template<class Thermo>
130 (
131  const scalar p,
132  const scalar T
133 ) const
134 {
135  scalar Cv_ = this->Cv(p, T);
136  return mu(p, T)*Cv_*(1.32 + 1.77*this->R()/Cv_);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
141 
142 template<class Thermo>
143 inline void Foam::sutherlandTransport<Thermo>::operator+=
144 (
146 )
147 {
148  scalar Y1 = this->Y();
149 
150  Thermo::operator+=(st);
151 
152  if (mag(this->Y()) > small)
153  {
154  Y1 /= this->Y();
155  scalar Y2 = st.Y()/this->Y();
156 
157  As_ = Y1*As_ + Y2*st.As_;
158  Ts_ = Y1*Ts_ + Y2*st.Ts_;
159  }
160 }
161 
162 
163 template<class Thermo>
164 inline void Foam::sutherlandTransport<Thermo>::operator*=
165 (
166  const scalar s
167 )
168 {
169  Thermo::operator*=(s);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
174 
175 template<class Thermo>
176 inline Foam::sutherlandTransport<Thermo> Foam::operator+
177 (
178  const sutherlandTransport<Thermo>& st1,
179  const sutherlandTransport<Thermo>& st2
180 )
181 {
182  Thermo t
183  (
184  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
185  );
186 
187  if (mag(t.Y()) < small)
188  {
190  (
191  t,
192  0,
193  st1.As_,
194  st1.Ts_
195  );
196  }
197  else
198  {
199  scalar Y1 = st1.Y()/t.Y();
200  scalar Y2 = st2.Y()/t.Y();
201 
203  (
204  t,
205  Y1*st1.As_ + Y2*st2.As_,
206  Y1*st1.Ts_ + Y2*st2.Ts_
207  );
208  }
209 }
210 
211 
212 template<class Thermo>
213 inline Foam::sutherlandTransport<Thermo> Foam::operator*
214 (
215  const scalar s,
217 )
218 {
220  (
221  s*static_cast<const Thermo&>(st),
222  st.As_,
223  st.Ts_
224  );
225 }
226 
227 
228 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
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))
A class for handling words, derived from string.
Definition: word.H:59
sutherlandTransport(const Thermo &t, const scalar As, const scalar Ts)
Construct from components.
autoPtr< sutherlandTransport > clone() const
Construct and return a clone.
static autoPtr< sutherlandTransport > New(const dictionary &dict)
const dimensionedScalar mu
Atomic mass unit.
const volScalarField & T
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
#define R(A, B, C, D, E, F, K, M)
PtrList< volScalarField > & Y
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
Transport package using Sutherland&#39;s formula.