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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "specie.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class Thermo>
34 (
35  const scalar mu1, const scalar T1,
36  const scalar mu2, const scalar T2
37 )
38 {
39  scalar rootT1 = sqrt(T1);
40  scalar mu1rootT2 = mu1*sqrt(T2);
41  scalar mu2rootT1 = mu2*rootT1;
42 
43  Ts_ = (mu2rootT1 - mu1rootT2)/(mu1rootT2/T1 - mu2rootT1/T2);
44 
45  As_ = mu1*(1.0 + Ts_/T1)/rootT1;
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class Thermo>
53 (
54  const Thermo& t,
55  const scalar As,
56  const scalar Ts
57 )
58 :
59  Thermo(t),
60  As_(As),
61  Ts_(Ts)
62 {}
63 
64 
65 template<class Thermo>
67 (
68  const Thermo& t,
69  const scalar mu1, const scalar T1,
70  const scalar mu2, const scalar T2
71 )
72 :
73  Thermo(t)
74 {
75  calcCoeffs(mu1, T1, mu2, T2);
76 }
77 
78 
79 template<class Thermo>
81 (
82  const word& name,
83  const sutherlandTransport& st
84 )
85 :
86  Thermo(name, st),
87  As_(st.As_),
88  Ts_(st.Ts_)
89 {}
90 
91 
92 template<class Thermo>
95 {
97 }
98 
99 
100 template<class Thermo>
103 (
104  const dictionary& dict
105 )
106 {
108 }
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
113 template<class Thermo>
114 inline Foam::scalar Foam::sutherlandTransport<Thermo>::mu
115 (
116  const scalar p,
117  const scalar T
118 ) const
119 {
120  return As_*::sqrt(T)/(1.0 + Ts_/T);
121 }
122 
123 
124 template<class Thermo>
126 (
127  const scalar p, const scalar T
128 ) const
129 {
130  scalar Cv_ = this->Cv(p, T);
131  return mu(p, T)*Cv_*(1.32 + 1.77*this->R()/Cv_);
132 }
133 
134 
135 template<class Thermo>
137 (
138  const scalar p,
139  const scalar T
140 ) const
141 {
142 
143  return kappa(p, T)/this->Cp(p, T);
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
148 
149 template<class Thermo>
151 (
153 )
154 {
155  scalar Y1 = this->Y();
156 
157  Thermo::operator+=(st);
158 
159  if (mag(this->Y()) > SMALL)
160  {
161  Y1 /= this->Y();
162  scalar Y2 = st.Y()/this->Y();
163 
164  As_ = Y1*As_ + Y2*st.As_;
165  Ts_ = Y1*Ts_ + Y2*st.Ts_;
166  }
167 }
168 
169 
170 template<class Thermo>
172 (
173  const scalar s
174 )
175 {
176  Thermo::operator*=(s);
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
181 
182 template<class Thermo>
183 inline Foam::sutherlandTransport<Thermo> Foam::operator+
184 (
185  const sutherlandTransport<Thermo>& st1,
186  const sutherlandTransport<Thermo>& st2
187 )
188 {
189  Thermo t
190  (
191  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
192  );
193 
194  if (mag(t.Y()) < SMALL)
195  {
197  (
198  t,
199  st1.As_,
200  st1.Ts_
201  );
202  }
203  else
204  {
205  scalar Y1 = st1.Y()/t.Y();
206  scalar Y2 = st2.Y()/t.Y();
207 
208  return sutherlandTransport<Thermo>
209  (
210  t,
211  Y1*st1.As_ + Y2*st2.As_,
212  Y1*st1.Ts_ + Y2*st2.Ts_
213  );
214  }
215 }
216 
217 
218 template<class Thermo>
219 inline Foam::sutherlandTransport<Thermo> Foam::operator*
220 (
221  const scalar s,
222  const sutherlandTransport<Thermo>& st
223 )
224 {
225  return sutherlandTransport<Thermo>
226  (
227  s*static_cast<const Thermo&>(st),
228  st.As_,
229  st.Ts_
230  );
231 }
232 
233 
234 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::sutherlandTransport::kappa
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
Definition: sutherlandTransportI.H:126
Foam::sutherlandTransport::alphah
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
Definition: sutherlandTransportI.H:137
Cv
const volScalarField & Cv
Definition: EEqn.H:8
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::mu
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
Definition: sutherlandTransportI.H:115
Foam::sutherlandTransport::clone
autoPtr< sutherlandTransport > clone() const
Construct and return a clone.
Definition: sutherlandTransportI.H:94
Foam::sutherlandTransport
Transport package using Sutherland's formula.
Definition: sutherlandTransport.H:59
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::sutherlandTransport::New
static autoPtr< sutherlandTransport > New(const dictionary &dict)
Definition: sutherlandTransportI.H:103
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::sutherlandTransport::sutherlandTransport
sutherlandTransport(const Thermo &t, const scalar As, const scalar Ts)
Construct from components.
Definition: sutherlandTransportI.H:53
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59