constTransportI.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-2015 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class Thermo>
30 (
31  const Thermo& t,
32  const scalar mu,
33  const scalar Pr
34 )
35 :
36  Thermo(t),
37  mu_(mu),
38  rPr_(1.0/Pr)
39 {}
40 
41 
42 template<class Thermo>
44 (
45  const word& name,
46  const constTransport& ct
47 )
48 :
49  Thermo(name, ct),
50  mu_(ct.mu_),
51  rPr_(ct.rPr_)
52 {}
53 
54 
55 template<class Thermo>
58 {
60  (
61  new constTransport<Thermo>(*this)
62  );
63 }
64 
65 
66 template<class Thermo>
69 (
70  Istream& is
71 )
72 {
74  (
76  );
77 }
78 
79 
80 template<class Thermo>
83 (
84  const dictionary& dict
85 )
86 {
88  (
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class Thermo>
97 inline Foam::scalar Foam::constTransport<Thermo>::mu
98 (
99  const scalar p,
100  const scalar T
101 ) const
102 {
103  return mu_;
104 }
105 
106 
107 template<class Thermo>
108 inline Foam::scalar Foam::constTransport<Thermo>::kappa
109 (
110  const scalar p,
111  const scalar T
112 ) const
113 {
114  return this->Cp(p, T)*mu(p, T)*rPr_;
115 }
116 
117 
118 template<class Thermo>
119 inline Foam::scalar Foam::constTransport<Thermo>::alphah
120 (
121  const scalar p,
122  const scalar T
123 ) const
124 {
125  return mu(p, T)*rPr_;
126 }
127 
128 
129 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
130 
131 template<class Thermo>
133 (
134  const constTransport<Thermo>& ct
135 )
136 {
137  Thermo::operator=(ct);
138 
139  mu_ = ct.mu_;
140  rPr_ = ct.rPr_;
141 
142  return *this;
143 }
144 
145 
146 template<class Thermo>
148 (
149  const constTransport<Thermo>& st
150 )
151 {
152  scalar molr1 = this->nMoles();
153 
154  Thermo::operator+=(st);
155 
156  if (mag(molr1) + mag(st.nMoles()) > SMALL)
157  {
158  molr1 /= this->nMoles();
159  scalar molr2 = st.nMoles()/this->nMoles();
160 
161  mu_ = molr1*mu_ + molr2*st.mu_;
162  rPr_ = 1.0/(molr1/rPr_ + molr2/st.rPr_);
163  }
164 }
165 
166 
167 template<class Thermo>
169 (
170  const constTransport<Thermo>& st
171 )
172 {
173  scalar molr1 = this->nMoles();
174 
175  Thermo::operator-=(st);
176 
177  if (mag(molr1) + mag(st.nMoles()) > SMALL)
178  {
179  molr1 /= this->nMoles();
180  scalar molr2 = st.nMoles()/this->nMoles();
181 
182  mu_ = molr1*mu_ - molr2*st.mu_;
183  rPr_ = 1.0/(molr1/rPr_ - molr2/st.rPr_);
184  }
185 }
186 
187 
188 template<class Thermo>
190 (
191  const scalar s
192 )
193 {
194  Thermo::operator*=(s);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
199 
200 template<class Thermo>
201 inline Foam::constTransport<Thermo> Foam::operator+
202 (
203  const constTransport<Thermo>& ct1,
204  const constTransport<Thermo>& ct2
205 )
206 {
207  Thermo t
208  (
209  static_cast<const Thermo&>(ct1) + static_cast<const Thermo&>(ct2)
210  );
211 
212  if (mag(ct1.nMoles()) + mag(ct2.nMoles()) < SMALL)
213  {
215  (
216  t,
217  0,
218  ct1.rPr_
219  );
220  }
221  else
222  {
223  scalar molr1 = ct1.nMoles()/t.nMoles();
224  scalar molr2 = ct2.nMoles()/t.nMoles();
225 
226  return constTransport<Thermo>
227  (
228  t,
229  molr1*ct1.mu_ + molr2*ct2.mu_,
230  1.0/(molr1/ct1.rPr_ + molr2/ct2.rPr_)
231  );
232  }
233 }
234 
235 
236 template<class Thermo>
237 inline Foam::constTransport<Thermo> Foam::operator-
238 (
239  const constTransport<Thermo>& ct1,
240  const constTransport<Thermo>& ct2
241 )
242 {
243  Thermo t
244  (
245  static_cast<const Thermo&>(ct1) - static_cast<const Thermo&>(ct2)
246  );
247 
248  if (mag(ct1.nMoles()) + mag(ct2.nMoles()) < SMALL)
249  {
250  return constTransport<Thermo>
251  (
252  t,
253  0,
254  ct1.rPr_
255  );
256  }
257  else
258  {
259  scalar molr1 = ct1.nMoles()/t.nMoles();
260  scalar molr2 = ct2.nMoles()/t.nMoles();
261 
262  return constTransport<Thermo>
263  (
264  t,
265  molr1*ct1.mu_ - molr2*ct2.mu_,
266  1.0/(molr1/ct1.rPr_ - molr2/ct2.rPr_)
267  );
268  }
269 }
270 
271 
272 template<class Thermo>
273 inline Foam::constTransport<Thermo> Foam::operator*
274 (
275  const scalar s,
276  const constTransport<Thermo>& ct
277 )
278 {
279  return constTransport<Thermo>
280  (
281  s*static_cast<const Thermo&>(ct),
282  ct.mu_,
283  1.0/ct.rPr_
284  );
285 }
286 
287 
288 template<class Thermo>
289 inline Foam::constTransport<Thermo> Foam::operator==
290 (
291  const constTransport<Thermo>& ct1,
292  const constTransport<Thermo>& ct2
293 )
294 {
295  return ct2 - ct1;
296 }
297 
298 
299 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constTransport::rPr_
scalar rPr_
Reciprocal Prandtl Number [].
Definition: constTransport.H:100
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::constTransport::New
static autoPtr< constTransport > New(Istream &is)
Definition: constTransportI.H:69
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::constTransport::mu_
scalar mu_
Constant dynamic viscosity [Pa.s].
Definition: constTransport.H:97
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::constTransport::kappa
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
Definition: constTransportI.H:109
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::constTransport
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
Definition: constTransport.H:47
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::constTransport::alphah
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
Definition: constTransportI.H:120
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::constTransport::constTransport
constTransport(const Thermo &t, const scalar mu, const scalar Pr)
Construct from components.
Definition: constTransportI.H:30
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::constTransport::mu
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
Definition: constTransportI.H:98
Foam::constTransport::clone
autoPtr< constTransport > clone() const
Construct and return a clone.
Definition: constTransportI.H:57
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47