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-2023 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 bool constPr,
34  const scalar rPr,
35  const scalar kappa
36 )
37 :
38  Thermo(t),
39  mu_(mu),
40  constPr_(constPr),
41  rPr_(rPr),
42  kappa_(kappa)
43 {}
44 
45 
46 template<class Thermo>
48 (
49  const word& name,
50  const constTransport& ct
51 )
52 :
53  Thermo(name, ct),
54  mu_(ct.mu_),
55  constPr_(ct.constPr_),
56  rPr_(ct.rPr_),
57  kappa_(ct.kappa_)
58 {}
59 
60 
61 template<class Thermo>
64 {
66  (
67  new constTransport<Thermo>(*this)
68  );
69 }
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class Thermo>
76 (
77  const scalar p,
78  const scalar T
79 ) const
80 {
81  return mu_;
82 }
83 
84 
85 template<class Thermo>
87 (
88  const scalar p,
89  const scalar T
90 ) const
91 {
92  return constPr_ ? this->Cp(p, T)*mu(p, T)*rPr_ : kappa_;
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
97 
98 template<class Thermo>
100 (
101  const constTransport<Thermo>& st
102 )
103 {
104  scalar Y1 = this->Y();
105 
106  Thermo::operator+=(st);
107 
108  if (mag(this->Y()) > small)
109  {
110  if
111  (
113  && (constPr_ != st.constPr_)
114  )
115  {
117  << "Constant " << (constPr_ ? "Pr" : "kappa") << " for "
118  << (this->name().size() ? this->name() : "others") << " but "
119  << "constant " << (st.constPr_ ? "Pr" : "kappa") << " for "
120  << (st.name().size() ? st.name() : "others")
121  << exit(FatalError);
122  }
123 
124  Y1 /= this->Y();
125  scalar Y2 = st.Y()/this->Y();
126 
127  mu_ = Y1*mu_ + Y2*st.mu_;
128  rPr_ = constPr_ ? 1/(Y1/rPr_ + Y2/st.rPr_) : NaN;
129  kappa_ = constPr_ ? NaN : Y1*kappa_ + Y2*st.kappa_;
130  }
131 }
132 
133 
134 template<class Thermo>
136 (
137  const scalar s
138 )
139 {
140  Thermo::operator*=(s);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
145 
146 template<class Thermo>
147 inline Foam::constTransport<Thermo> Foam::operator+
148 (
149  const constTransport<Thermo>& ct1,
150  const constTransport<Thermo>& ct2
151 )
152 {
153  Thermo t
154  (
155  static_cast<const Thermo&>(ct1) + static_cast<const Thermo&>(ct2)
156  );
157 
158  if (mag(t.Y()) < small)
159  {
161  (
162  t,
163  0,
164  ct1.constPr_,
165  ct1.rPr_,
166  0
167  );
168  }
169  else
170  {
171  scalar Y1 = ct1.Y()/t.Y();
172  scalar Y2 = ct2.Y()/t.Y();
173 
174  if
175  (
176  constTransport<Thermo>::debug
177  && (ct1.constPr_ != ct2.constPr_)
178  )
179  {
181  << "Constant " << (ct1.constPr_ ? "Pr" : "kappa") << " for "
182  << (ct1.name().size() ? ct1.name() : "others") << " but "
183  << "constant " << (ct2.constPr_ ? "Pr" : "kappa") << " for "
184  << (ct2.name().size() ? ct2.name() : "others")
185  << exit(FatalError);
186  }
187 
188  return constTransport<Thermo>
189  (
190  t,
191  Y1*ct1.mu_ + Y2*ct2.mu_,
192  ct1.constPr_,
193  ct1.constPr_ ? 1/(Y1/ct1.rPr_ + Y2/ct2.rPr_) : NaN,
194  ct1.constPr_ ? NaN : Y1*ct1.kappa_ + Y2*ct2.kappa_
195  );
196  }
197 }
198 
199 
200 template<class Thermo>
201 inline Foam::constTransport<Thermo> Foam::operator*
202 (
203  const scalar s,
204  const constTransport<Thermo>& ct
205 )
206 {
207  return constTransport<Thermo>
208  (
209  s*static_cast<const Thermo&>(ct),
210  ct.mu_,
211  ct.constPr_,
212  ct.rPr_,
213  ct.kappa_
214  );
215 }
216 
217 
218 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Transport package with constant properties.
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
scalar kappa(const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
autoPtr< constTransport > clone() const
Construct and return a clone.
constTransport(const Thermo &t, const scalar mu, const bool constPr, const scalar rPr, const scalar kappa)
Construct from components.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
static Type NaN()
Return a primitive with all components set to NaN.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p
PtrList< volScalarField > & Y