All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
powerLawTransportI.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 
13 License
14  This file is part of QGDsolver, based on OpenFOAM library.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Class
30  Foam::powerLawTransport
31 
32 Group
33  grpPowerLawQGD
34 
35 Description:
36  Definition of templates
37 
38 SourceFiles:
39  powerLawTransportI.H
40  powerLawTransport.C
41 
42 \*----------------------------------------------------------------------------*/
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Thermo>
48 (
49  const Thermo& t,
50  const scalar mu0,
51  const scalar T0,
52  const scalar kk,
53  const scalar Pr
54 )
55 :
56  Thermo(t),
57  mu0_(mu0),
58  T0_(T0),
59  k_(kk),
60  rPr_(1.0/Pr)
61 {}
62 
63 
64 template<class Thermo>
66 (
67  const word& name,
68  const powerLawTransport& ct
69 )
70 :
71  Thermo(name, ct),
72  mu0_(ct.mu0_),
73  T0_(ct.T0_),
74  k_(ct.k_),
75  rPr_(ct.rPr_)
76 {}
77 
78 
79 template<class Thermo>
80 inline Foam::autoPtr<Foam::powerLawTransport<Thermo> >
82 {
83  return autoPtr<powerLawTransport<Thermo> >
84  (
86  );
87 }
88 
89 
90 template<class Thermo>
91 inline Foam::autoPtr<Foam::powerLawTransport<Thermo> >
93 (
94  Istream& is
95 )
96 {
97  return autoPtr<powerLawTransport<Thermo> >
98  (
100  );
101 }
102 
103 
104 template<class Thermo>
105 inline Foam::autoPtr<Foam::powerLawTransport<Thermo> >
107 (
108  const dictionary& dict
109 )
110 {
111  return autoPtr<powerLawTransport<Thermo> >
112  (
113  new powerLawTransport<Thermo>(dict)
114  );
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class Thermo>
121 inline Foam::scalar Foam::powerLawTransport<Thermo>::mu
122 (
123  const scalar p,
124  const scalar T
125 ) const
126 {
127  return mu0_ * pow(T / T0_, k_);
128 }
129 
130 
131 template<class Thermo>
132 inline Foam::scalar Foam::powerLawTransport<Thermo>::kappa
133 (
134  const scalar p,
135  const scalar T
136 ) const
137 {
138  return this->Cp(p, T)*mu(p, T)*rPr_;
139 }
140 
141 
142 template<class Thermo>
143 inline Foam::scalar Foam::powerLawTransport<Thermo>::alphah
144 (
145  const scalar p,
146  const scalar T
147 ) const
148 {
149  return mu(p, T)*rPr_;
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
155 template<class Thermo>
157 (
158  const powerLawTransport<Thermo>& ct
159 )
160 {
161  Thermo::operator=(ct);
162 
163  mu0_ = ct.mu0_;
164  T0_ = ct.T0_;
165  k_ = ct.k_;
166  rPr_ = ct.rPr_;
167 
168  return *this;
169 }
170 
171 
172 template<class Thermo>
173 inline void Foam::powerLawTransport<Thermo>::operator+=
174 (
175  const powerLawTransport<Thermo>& st
176 )
177 {
178  scalar molr1 = this->nMoles();
179 
180  Thermo::operator+=(st);
181 
182  if (mag(molr1) + mag(st.nMoles()) > SMALL)
183  {
184  molr1 /= this->nMoles();
185  scalar molr2 = st.nMoles()/this->nMoles();
187  rPr_ = 1.0/(molr1/rPr_ + molr2/st.rPr_);
188  }
189 }
190 
191 
192 template<class Thermo>
193 inline void Foam::powerLawTransport<Thermo>::operator-=
194 (
195  const powerLawTransport<Thermo>& st
196 )
197 {
198  scalar molr1 = this->nMoles();
199 
200  Thermo::operator-=(st);
201 
202  if (mag(molr1) + mag(st.nMoles()) > SMALL)
203  {
204  molr1 /= this->nMoles();
205  scalar molr2 = st.nMoles()/this->nMoles();
207  rPr_ = 1.0/(molr1/rPr_ - molr2/st.rPr_);
208  }
209 }
210 
211 
212 template<class Thermo>
213 inline void Foam::powerLawTransport<Thermo>::operator*=
214 (
215  const scalar s
216 )
217 {
218  Thermo::operator*=(s);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
223 
224 template<class Thermo>
225 inline Foam::powerLawTransport<Thermo> Foam::operator+
226 (
227  const powerLawTransport<Thermo>& ct1,
228  const powerLawTransport<Thermo>& ct2
229 )
230 {
231  notImplemented ("inline Foam::powerLawTransport<Thermo> Foam::operator+");
232 }
233 
234 
235 template<class Thermo>
236 inline Foam::powerLawTransport<Thermo> Foam::operator-
237 (
238  const powerLawTransport<Thermo>& ct1,
239  const powerLawTransport<Thermo>& ct2
240 )
241 {
242  notImplemented ("inline Foam::powerLawTransport<Thermo> Foam::operator-");
243 }
244 
245 
246 template<class Thermo>
247 inline Foam::powerLawTransport<Thermo> Foam::operator*
248 (
249  const scalar s,
250  const powerLawTransport<Thermo>& ct
251 )
252 {
253  notImplemented ("inline Foam::powerLawTransport<Thermo> Foam::operator*");
254 }
255 
256 
257 template<class Thermo>
258 inline Foam::powerLawTransport<Thermo> Foam::operator==
259 (
260  const powerLawTransport<Thermo>& ct1,
261  const powerLawTransport<Thermo>& ct2
262 )
263 {
264  notImplemented ("inline Foam::powerLawTransport<Thermo> Foam::operator==");
265 }
266 
267 
268 // ************************************************************************* //
volScalarField T0("T0", T)
autoPtr< powerLawTransport > clone() const
Construct and return a clone.
static autoPtr< powerLawTransport > New(Istream &is)
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/mK].
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
volScalarField & T
Definition: createFields.H:53
volScalarField & p
Definition: createFields.H:52
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...