constTransportI.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 // * * * * * * * * * * * * * * * * 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  const dictionary& dict
71 )
72 {
74  (
76  );
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class Thermo>
83 inline Foam::scalar Foam::constTransport<Thermo>::mu
84 (
85  const scalar p,
86  const scalar T
87 ) const
88 {
89  return mu_;
90 }
91 
92 
93 template<class Thermo>
94 inline Foam::scalar Foam::constTransport<Thermo>::kappa
95 (
96  const scalar p,
97  const scalar T
98 ) const
99 {
100  return this->Cp(p, T)*mu(p, T)*rPr_;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
105 
106 template<class Thermo>
107 inline void Foam::constTransport<Thermo>::operator+=
108 (
109  const constTransport<Thermo>& st
110 )
111 {
112  scalar Y1 = this->Y();
113 
114  Thermo::operator+=(st);
115 
116  if (mag(this->Y()) > small)
117  {
118  Y1 /= this->Y();
119  scalar Y2 = st.Y()/this->Y();
120 
121  mu_ = Y1*mu_ + Y2*st.mu_;
122  rPr_ = 1.0/(Y1/rPr_ + Y2/st.rPr_);
123  }
124 }
125 
126 
127 template<class Thermo>
128 inline void Foam::constTransport<Thermo>::operator*=
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::constTransport<Thermo> Foam::operator+
141 (
142  const constTransport<Thermo>& ct1,
143  const constTransport<Thermo>& ct2
144 )
145 {
146  Thermo t
147  (
148  static_cast<const Thermo&>(ct1) + static_cast<const Thermo&>(ct2)
149  );
150 
151  if (mag(t.Y()) < small)
152  {
154  (
155  t,
156  0,
157  ct1.rPr_
158  );
159  }
160  else
161  {
162  scalar Y1 = ct1.Y()/t.Y();
163  scalar Y2 = ct2.Y()/t.Y();
164 
166  (
167  t,
168  Y1*ct1.mu_ + Y2*ct2.mu_,
169  1.0/(Y1/ct1.rPr_ + Y2/ct2.rPr_)
170  );
171  }
172 }
173 
174 
175 template<class Thermo>
176 inline Foam::constTransport<Thermo> Foam::operator*
177 (
178  const scalar s,
179  const constTransport<Thermo>& ct
180 )
181 {
183  (
184  s*static_cast<const Thermo&>(ct),
185  ct.mu_,
186  1.0/ct.rPr_
187  );
188 }
189 
190 
191 // ************************************************************************* //
dictionary dict
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
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
autoPtr< constTransport > clone() const
Construct and return a clone.
const dimensionedScalar mu
Atomic mass unit.
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
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
static autoPtr< constTransport > New(const dictionary &dict)
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].