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-2016 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>
132 inline void Foam::constTransport<Thermo>::operator=
133 (
134  const constTransport<Thermo>& ct
135 )
136 {
137  Thermo::operator=(ct);
138 
139  mu_ = ct.mu_;
140  rPr_ = ct.rPr_;
141 }
142 
143 
144 template<class Thermo>
145 inline void Foam::constTransport<Thermo>::operator+=
146 (
147  const constTransport<Thermo>& st
148 )
149 {
150  scalar molr1 = this->nMoles();
151 
152  Thermo::operator+=(st);
153 
154  if (mag(molr1) + mag(st.nMoles()) > SMALL)
155  {
156  molr1 /= this->nMoles();
157  scalar molr2 = st.nMoles()/this->nMoles();
158 
159  mu_ = molr1*mu_ + molr2*st.mu_;
160  rPr_ = 1.0/(molr1/rPr_ + molr2/st.rPr_);
161  }
162 }
163 
164 
165 template<class Thermo>
166 inline void Foam::constTransport<Thermo>::operator-=
167 (
168  const constTransport<Thermo>& st
169 )
170 {
171  scalar molr1 = this->nMoles();
172 
173  Thermo::operator-=(st);
174 
175  if (mag(molr1) + mag(st.nMoles()) > SMALL)
176  {
177  molr1 /= this->nMoles();
178  scalar molr2 = st.nMoles()/this->nMoles();
179 
180  mu_ = molr1*mu_ - molr2*st.mu_;
181  rPr_ = 1.0/(molr1/rPr_ - molr2/st.rPr_);
182  }
183 }
184 
185 
186 template<class Thermo>
187 inline void Foam::constTransport<Thermo>::operator*=
188 (
189  const scalar s
190 )
191 {
192  Thermo::operator*=(s);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
197 
198 template<class Thermo>
199 inline Foam::constTransport<Thermo> Foam::operator+
200 (
201  const constTransport<Thermo>& ct1,
202  const constTransport<Thermo>& ct2
203 )
204 {
205  Thermo t
206  (
207  static_cast<const Thermo&>(ct1) + static_cast<const Thermo&>(ct2)
208  );
209 
210  if (mag(ct1.nMoles()) + mag(ct2.nMoles()) < SMALL)
211  {
213  (
214  t,
215  0,
216  ct1.rPr_
217  );
218  }
219  else
220  {
221  scalar molr1 = ct1.nMoles()/t.nMoles();
222  scalar molr2 = ct2.nMoles()/t.nMoles();
223 
225  (
226  t,
227  molr1*ct1.mu_ + molr2*ct2.mu_,
228  1.0/(molr1/ct1.rPr_ + molr2/ct2.rPr_)
229  );
230  }
231 }
232 
233 
234 template<class Thermo>
235 inline Foam::constTransport<Thermo> Foam::operator-
236 (
237  const constTransport<Thermo>& ct1,
238  const constTransport<Thermo>& ct2
239 )
240 {
241  Thermo t
242  (
243  static_cast<const Thermo&>(ct1) - static_cast<const Thermo&>(ct2)
244  );
245 
246  if (mag(ct1.nMoles()) + mag(ct2.nMoles()) < SMALL)
247  {
249  (
250  t,
251  0,
252  ct1.rPr_
253  );
254  }
255  else
256  {
257  scalar molr1 = ct1.nMoles()/t.nMoles();
258  scalar molr2 = ct2.nMoles()/t.nMoles();
259 
261  (
262  t,
263  molr1*ct1.mu_ - molr2*ct2.mu_,
264  1.0/(molr1/ct1.rPr_ - molr2/ct2.rPr_)
265  );
266  }
267 }
268 
269 
270 template<class Thermo>
271 inline Foam::constTransport<Thermo> Foam::operator*
272 (
273  const scalar s,
274  const constTransport<Thermo>& ct
275 )
276 {
278  (
279  s*static_cast<const Thermo&>(ct),
280  ct.mu_,
281  1.0/ct.rPr_
282  );
283 }
284 
285 
286 template<class Thermo>
287 inline Foam::constTransport<Thermo> Foam::operator==
288 (
289  const constTransport<Thermo>& ct1,
290  const constTransport<Thermo>& ct2
291 )
292 {
293  return ct2 - ct1;
294 }
295 
296 
297 // ************************************************************************* //
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
autoPtr< constTransport > clone() const
Construct and return a clone.
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
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
static autoPtr< constTransport > New(Istream &is)
const dimensionedScalar mu
Atomic mass unit.
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:53