sutherlandTransportI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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  Istream& is
106 )
107 {
109  (
111  );
112 }
113 
114 
115 template<class Thermo>
118 (
119  const dictionary& dict
120 )
121 {
123  (
125  );
126 }
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class Thermo>
132 inline Foam::scalar Foam::sutherlandTransport<Thermo>::mu
133 (
134  const scalar p,
135  const scalar T
136 ) const
137 {
138  return As_*::sqrt(T)/(1.0 + Ts_/T);
139 }
140 
141 
142 template<class Thermo>
144 (
145  const scalar p, const scalar T
146 ) const
147 {
148  scalar Cv_ = this->Cv(p, T);
149  return mu(p, T)*Cv_*(1.32 + 1.77*this->R()/Cv_);
150 }
151 
152 
153 template<class Thermo>
155 (
156  const scalar p,
157  const scalar T
158 ) const
159 {
160 
161  return kappa(p, T)/this->Cpv(p, T);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
166 
167 template<class Thermo>
170 (
172 )
173 {
174  Thermo::operator=(st);
175 
176  As_ = st.As_;
177  Ts_ = st.Ts_;
178 
179  return *this;
180 }
181 
182 
183 template<class Thermo>
185 (
187 )
188 {
189  scalar molr1 = this->nMoles();
190 
191  Thermo::operator+=(st);
192 
193  molr1 /= this->nMoles();
194  scalar molr2 = st.nMoles()/this->nMoles();
195 
196  As_ = molr1*As_ + molr2*st.As_;
197  Ts_ = molr1*Ts_ + molr2*st.Ts_;
198 }
199 
200 
201 template<class Thermo>
203 (
205 )
206 {
207  scalar molr1 = this->nMoles();
208 
209  Thermo::operator-=(st);
210 
211  molr1 /= this->nMoles();
212  scalar molr2 = st.nMoles()/this->nMoles();
213 
214  As_ = molr1*As_ - molr2*st.As_;
215  Ts_ = molr1*Ts_ - molr2*st.Ts_;
216 }
217 
218 
219 template<class Thermo>
221 (
222  const scalar s
223 )
224 {
225  Thermo::operator*=(s);
226 }
227 
228 
229 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
230 
231 template<class Thermo>
232 inline Foam::sutherlandTransport<Thermo> Foam::operator+
233 (
234  const sutherlandTransport<Thermo>& st1,
235  const sutherlandTransport<Thermo>& st2
236 )
237 {
238  Thermo t
239  (
240  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
241  );
242 
243  scalar molr1 = st1.nMoles()/t.nMoles();
244  scalar molr2 = st2.nMoles()/t.nMoles();
245 
247  (
248  t,
249  molr1*st1.As_ + molr2*st2.As_,
250  molr1*st1.Ts_ + molr2*st2.Ts_
251  );
252 }
253 
254 
255 template<class Thermo>
256 inline Foam::sutherlandTransport<Thermo> Foam::operator-
257 (
258  const sutherlandTransport<Thermo>& st1,
259  const sutherlandTransport<Thermo>& st2
260 )
261 {
262  Thermo t
263  (
264  static_cast<const Thermo&>(st1) - static_cast<const Thermo&>(st2)
265  );
266 
267  scalar molr1 = st1.nMoles()/t.nMoles();
268  scalar molr2 = st2.nMoles()/t.nMoles();
269 
270  return sutherlandTransport<Thermo>
271  (
272  t,
273  molr1*st1.As_ - molr2*st2.As_,
274  molr1*st1.Ts_ - molr2*st2.Ts_
275  );
276 }
277 
278 
279 template<class Thermo>
280 inline Foam::sutherlandTransport<Thermo> Foam::operator*
281 (
282  const scalar s,
283  const sutherlandTransport<Thermo>& st
284 )
285 {
286  return sutherlandTransport<Thermo>
287  (
288  s*static_cast<const Thermo&>(st),
289  st.As_,
290  st.Ts_
291  );
292 }
293 
294 
295 template<class Thermo>
296 inline Foam::sutherlandTransport<Thermo> Foam::operator==
297 (
298  const sutherlandTransport<Thermo>& st1,
299  const sutherlandTransport<Thermo>& st2
300 )
301 {
302  return st2 - st1;
303 }
304 
305 
306 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::sutherlandTransport::kappa
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
Definition: sutherlandTransportI.H:144
Foam::sutherlandTransport::alphah
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
Definition: sutherlandTransportI.H:155
Foam::sutherlandTransport::Ts_
scalar Ts_
Definition: sutherlandTransport.H:103
specie.H
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
R
#define R(A, B, C, D, E, F, K, M)
Foam::sutherlandTransport::As_
scalar As_
Definition: sutherlandTransport.H:103
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::sutherlandTransport::mu
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
Definition: sutherlandTransportI.H:133
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::sutherlandTransport::clone
autoPtr< sutherlandTransport > clone() const
Construct and return a clone.
Definition: sutherlandTransportI.H:92
Foam::sutherlandTransport
Transport package using Sutherland's formula.
Definition: sutherlandTransport.H:53
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
s
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))
Foam::sutherlandTransport::calcCoeffs
void calcCoeffs(const scalar mu1, const scalar T1, const scalar mu2, const scalar T2)
Calculate the Sutherland coefficients.
Definition: sutherlandTransportI.H:32
Foam::sutherlandTransport::New
static autoPtr< sutherlandTransport > New(Istream &is)
Definition: sutherlandTransportI.H:104
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::sutherlandTransport::sutherlandTransport
sutherlandTransport(const Thermo &t, const scalar As, const scalar Ts)
Construct from components.
Definition: sutherlandTransportI.H:51
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47