sutherlandTransportI.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-2019 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  const dictionary& dict
106 )
107 {
109  (
111  );
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class Thermo>
118 inline Foam::scalar Foam::sutherlandTransport<Thermo>::mu
119 (
120  const scalar p,
121  const scalar T
122 ) const
123 {
124  return As_*::sqrt(T)/(1.0 + Ts_/T);
125 }
126 
127 
128 template<class Thermo>
130 (
131  const scalar p, const scalar T
132 ) const
133 {
134  scalar Cv_ = this->Cv(p, T);
135  return mu(p, T)*Cv_*(1.32 + 1.77*this->R()/Cv_);
136 }
137 
138 
139 template<class Thermo>
141 (
142  const scalar p,
143  const scalar T
144 ) const
145 {
146 
147  return kappa(p, T)/this->Cp(p, T);
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
152 
153 template<class Thermo>
154 inline void Foam::sutherlandTransport<Thermo>::operator+=
155 (
157 )
158 {
159  scalar Y1 = this->Y();
160 
161  Thermo::operator+=(st);
162 
163  if (mag(this->Y()) > small)
164  {
165  Y1 /= this->Y();
166  scalar Y2 = st.Y()/this->Y();
167 
168  As_ = Y1*As_ + Y2*st.As_;
169  Ts_ = Y1*Ts_ + Y2*st.Ts_;
170  }
171 }
172 
173 
174 template<class Thermo>
175 inline void Foam::sutherlandTransport<Thermo>::operator*=
176 (
177  const scalar s
178 )
179 {
180  Thermo::operator*=(s);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
185 
186 template<class Thermo>
187 inline Foam::sutherlandTransport<Thermo> Foam::operator+
188 (
189  const sutherlandTransport<Thermo>& st1,
190  const sutherlandTransport<Thermo>& st2
191 )
192 {
193  Thermo t
194  (
195  static_cast<const Thermo&>(st1) + static_cast<const Thermo&>(st2)
196  );
197 
198  if (mag(t.Y()) < small)
199  {
201  (
202  t,
203  0,
204  st1.As_,
205  st1.Ts_
206  );
207  }
208  else
209  {
210  scalar Y1 = st1.Y()/t.Y();
211  scalar Y2 = st2.Y()/t.Y();
212 
214  (
215  t,
216  Y1*st1.As_ + Y2*st2.As_,
217  Y1*st1.Ts_ + Y2*st2.Ts_
218  );
219  }
220 }
221 
222 
223 template<class Thermo>
224 inline Foam::sutherlandTransport<Thermo> Foam::operator*
225 (
226  const scalar s,
228 )
229 {
231  (
232  s*static_cast<const Thermo&>(st),
233  st.As_,
234  st.Ts_
235  );
236 }
237 
238 
239 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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
sutherlandTransport(const Thermo &t, const scalar As, const scalar Ts)
Construct from components.
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/m/s].
autoPtr< sutherlandTransport > clone() const
Construct and return a clone.
static autoPtr< sutherlandTransport > New(const dictionary &dict)
const volScalarField & T
const dimensionedScalar mu
Atomic mass unit.
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
#define R(A, B, C, D, E, F, K, M)
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
Transport package using Sutherland&#39;s formula.