thermophysicalModels/specie/transport: Added WLF transport model
Description
Transport package using the Williams-Landel-Ferry model.
Templated into a given thermodynamics package (needed for thermal
conductivity).
Dynamic viscosity [kg/m.s]
\f[
\mu = \mu_0 \exp \left(\frac{-C_1 ( T - T_r )}{C_2 + T - T_r}\right)
\f]
References:
\verbatim
Williams, M. L., Landel, R. F., & Ferry, J. D. (1955).
The temperature dependence of relaxation mechanisms
in amorphous polymers and other glass-forming liquids.
Journal of the American Chemical society, 77(14), 3701-3707.
\endverbatim
This commit is contained in:
@ -43,6 +43,7 @@ License
|
|||||||
|
|
||||||
#include "constTransport.H"
|
#include "constTransport.H"
|
||||||
#include "sutherlandTransport.H"
|
#include "sutherlandTransport.H"
|
||||||
|
#include "WLFTransport.H"
|
||||||
|
|
||||||
#include "icoPolynomial.H"
|
#include "icoPolynomial.H"
|
||||||
#include "hPolynomialThermo.H"
|
#include "hPolynomialThermo.H"
|
||||||
@ -408,6 +409,18 @@ makeThermos
|
|||||||
specie
|
specie
|
||||||
);
|
);
|
||||||
|
|
||||||
|
makeThermos
|
||||||
|
(
|
||||||
|
rhoThermo,
|
||||||
|
heRhoThermo,
|
||||||
|
pureMixture,
|
||||||
|
WLFTransport,
|
||||||
|
sensibleInternalEnergy,
|
||||||
|
hConstThermo,
|
||||||
|
rhoConst,
|
||||||
|
specie
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
|||||||
92
src/thermophysicalModels/specie/transport/WLF/WLFTransport.C
Normal file
92
src/thermophysicalModels/specie/transport/WLF/WLFTransport.C
Normal file
@ -0,0 +1,92 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration | Website: https://openfoam.org
|
||||||
|
\\ / A nd | Copyright (C) 2018 OpenFOAM Foundation
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
License
|
||||||
|
This file is part of OpenFOAM.
|
||||||
|
|
||||||
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
||||||
|
under the terms of the GNU General Public License as published by
|
||||||
|
the Free Software Foundation, either version 3 of the License, or
|
||||||
|
(at your option) any later version.
|
||||||
|
|
||||||
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
||||||
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||||
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
||||||
|
for more details.
|
||||||
|
|
||||||
|
You should have received a copy of the GNU General Public License
|
||||||
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#include "WLFTransport.H"
|
||||||
|
#include "IOstreams.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
Foam::scalar Foam::WLFTransport<Thermo>::readCoeff
|
||||||
|
(
|
||||||
|
const word& coeffName,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
{
|
||||||
|
return readScalar(dict.subDict("transport").lookup(coeffName));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
Foam::WLFTransport<Thermo>::WLFTransport(const dictionary& dict)
|
||||||
|
:
|
||||||
|
Thermo(dict),
|
||||||
|
mu0_(readCoeff("mu0", dict)),
|
||||||
|
Tr_(readCoeff("Tr", dict)),
|
||||||
|
C1_(readCoeff("C1", dict)),
|
||||||
|
C2_(readCoeff("C2", dict)),
|
||||||
|
rPr_(1.0/readScalar(dict.subDict("transport").lookup("Pr")))
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
void Foam::WLFTransport<Thermo>::write(Ostream& os) const
|
||||||
|
{
|
||||||
|
os << this->specie::name() << endl
|
||||||
|
<< token::BEGIN_BLOCK << incrIndent << nl;
|
||||||
|
|
||||||
|
Thermo::write(os);
|
||||||
|
|
||||||
|
dictionary dict("transport");
|
||||||
|
dict.add("mu0", mu0_);
|
||||||
|
dict.add("Tr", Tr_);
|
||||||
|
dict.add("C1", C1_);
|
||||||
|
dict.add("C2", C2_);
|
||||||
|
dict.add("Pr", 1.0/rPr_);
|
||||||
|
|
||||||
|
os << indent << dict.dictName() << dict
|
||||||
|
<< decrIndent << token::END_BLOCK << nl;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
Foam::Ostream& Foam::operator<<
|
||||||
|
(
|
||||||
|
Ostream& os,
|
||||||
|
const WLFTransport<Thermo>& wlft
|
||||||
|
)
|
||||||
|
{
|
||||||
|
wlft.write(os);
|
||||||
|
return os;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
210
src/thermophysicalModels/specie/transport/WLF/WLFTransport.H
Normal file
210
src/thermophysicalModels/specie/transport/WLF/WLFTransport.H
Normal file
@ -0,0 +1,210 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration | Website: https://openfoam.org
|
||||||
|
\\ / A nd | Copyright (C) 2018 OpenFOAM Foundation
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
License
|
||||||
|
This file is part of OpenFOAM.
|
||||||
|
|
||||||
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
||||||
|
under the terms of the GNU General Public License as published by
|
||||||
|
the Free Software Foundation, either version 3 of the License, or
|
||||||
|
(at your option) any later version.
|
||||||
|
|
||||||
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
||||||
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||||
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
||||||
|
for more details.
|
||||||
|
|
||||||
|
You should have received a copy of the GNU General Public License
|
||||||
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
|
||||||
|
Class
|
||||||
|
Foam::WLFTransport
|
||||||
|
|
||||||
|
Description
|
||||||
|
Transport package using the Williams-Landel-Ferry model.
|
||||||
|
|
||||||
|
Templated into a given thermodynamics package (needed for thermal
|
||||||
|
conductivity).
|
||||||
|
|
||||||
|
Dynamic viscosity [kg/m.s]
|
||||||
|
\f[
|
||||||
|
\mu = \mu_0 \exp \left(\frac{-C_1 ( T - T_r )}{C_2 + T - T_r}\right)
|
||||||
|
\f]
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
Williams, M. L., Landel, R. F., & Ferry, J. D. (1955).
|
||||||
|
The temperature dependence of relaxation mechanisms
|
||||||
|
in amorphous polymers and other glass-forming liquids.
|
||||||
|
Journal of the American Chemical society, 77(14), 3701-3707.
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
WLFTransportI.H
|
||||||
|
WLFTransport.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef WLFTransport_H
|
||||||
|
#define WLFTransport_H
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
// Forward declaration of friend functions and operators
|
||||||
|
|
||||||
|
template<class Thermo> class WLFTransport;
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline WLFTransport<Thermo> operator+
|
||||||
|
(
|
||||||
|
const WLFTransport<Thermo>&,
|
||||||
|
const WLFTransport<Thermo>&
|
||||||
|
);
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline WLFTransport<Thermo> operator*
|
||||||
|
(
|
||||||
|
const scalar,
|
||||||
|
const WLFTransport<Thermo>&
|
||||||
|
);
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
Ostream& operator<<
|
||||||
|
(
|
||||||
|
Ostream&,
|
||||||
|
const WLFTransport<Thermo>&
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class WLFTransport Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
class WLFTransport
|
||||||
|
:
|
||||||
|
public Thermo
|
||||||
|
{
|
||||||
|
// Private data
|
||||||
|
|
||||||
|
//- Dynamic viscosity at the reference temperature [Pa.s]
|
||||||
|
scalar mu0_;
|
||||||
|
|
||||||
|
//- Reference temperature [T]
|
||||||
|
scalar Tr_;
|
||||||
|
|
||||||
|
//- WLF coefficient 1 []
|
||||||
|
scalar C1_;
|
||||||
|
|
||||||
|
//- WLF coefficient 2 [T]
|
||||||
|
scalar C2_;
|
||||||
|
|
||||||
|
//- Reciprocal Prandtl Number []
|
||||||
|
scalar rPr_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Read coefficient from dictionary
|
||||||
|
scalar readCoeff(const word& coeffName, const dictionary& dict);
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct as named copy
|
||||||
|
inline WLFTransport(const word&, const WLFTransport&);
|
||||||
|
|
||||||
|
//- Construct from dictionary
|
||||||
|
WLFTransport(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Construct and return a clone
|
||||||
|
inline autoPtr<WLFTransport> clone() const;
|
||||||
|
|
||||||
|
// Selector from dictionary
|
||||||
|
inline static autoPtr<WLFTransport> New(const dictionary& dict);
|
||||||
|
|
||||||
|
|
||||||
|
// Member functions
|
||||||
|
|
||||||
|
//- Return the instantiated type name
|
||||||
|
static word typeName()
|
||||||
|
{
|
||||||
|
return "WLF<" + Thermo::typeName() + '>';
|
||||||
|
}
|
||||||
|
|
||||||
|
//- Dynamic viscosity [kg/ms]
|
||||||
|
inline scalar mu(const scalar p, const scalar T) const;
|
||||||
|
|
||||||
|
//- Thermal conductivity [W/mK]
|
||||||
|
inline scalar kappa(const scalar p, const scalar T) const;
|
||||||
|
|
||||||
|
//- Thermal diffusivity of enthalpy [kg/ms]
|
||||||
|
inline scalar alphah(const scalar p, const scalar T) const;
|
||||||
|
|
||||||
|
// Species diffusivity
|
||||||
|
// inline scalar D(const scalar p, const scalar T) const;
|
||||||
|
|
||||||
|
//- Write to Ostream
|
||||||
|
void write(Ostream& os) const;
|
||||||
|
|
||||||
|
|
||||||
|
// Member operators
|
||||||
|
|
||||||
|
inline void operator=(const WLFTransport&);
|
||||||
|
|
||||||
|
inline void operator+=(const WLFTransport&);
|
||||||
|
|
||||||
|
inline void operator*=(const scalar);
|
||||||
|
|
||||||
|
|
||||||
|
// Friend operators
|
||||||
|
|
||||||
|
friend WLFTransport operator+ <Thermo>
|
||||||
|
(
|
||||||
|
const WLFTransport&,
|
||||||
|
const WLFTransport&
|
||||||
|
);
|
||||||
|
|
||||||
|
friend WLFTransport operator* <Thermo>
|
||||||
|
(
|
||||||
|
const scalar,
|
||||||
|
const WLFTransport&
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Ostream Operator
|
||||||
|
|
||||||
|
friend Ostream& operator<< <Thermo>
|
||||||
|
(
|
||||||
|
Ostream&,
|
||||||
|
const WLFTransport&
|
||||||
|
);
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#include "WLFTransportI.H"
|
||||||
|
|
||||||
|
#ifdef NoRepository
|
||||||
|
#include "WLFTransport.C"
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
222
src/thermophysicalModels/specie/transport/WLF/WLFTransportI.H
Normal file
222
src/thermophysicalModels/specie/transport/WLF/WLFTransportI.H
Normal file
@ -0,0 +1,222 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration | Website: https://openfoam.org
|
||||||
|
\\ / A nd | Copyright (C) 2018 OpenFOAM Foundation
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
License
|
||||||
|
This file is part of OpenFOAM.
|
||||||
|
|
||||||
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
||||||
|
under the terms of the GNU General Public License as published by
|
||||||
|
the Free Software Foundation, either version 3 of the License, or
|
||||||
|
(at your option) any later version.
|
||||||
|
|
||||||
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
||||||
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||||
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
||||||
|
for more details.
|
||||||
|
|
||||||
|
You should have received a copy of the GNU General Public License
|
||||||
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#include "specie.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::WLFTransport<Thermo>::WLFTransport
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const WLFTransport& wlft
|
||||||
|
)
|
||||||
|
:
|
||||||
|
Thermo(name, wlft),
|
||||||
|
mu0_(wlft.mu0_),
|
||||||
|
Tr_(wlft.Tr_),
|
||||||
|
C1_(wlft.C1_),
|
||||||
|
C2_(wlft.C2_),
|
||||||
|
rPr_(wlft.rPr_)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::autoPtr<Foam::WLFTransport<Thermo>>
|
||||||
|
Foam::WLFTransport<Thermo>::clone() const
|
||||||
|
{
|
||||||
|
return autoPtr<WLFTransport<Thermo>>
|
||||||
|
(
|
||||||
|
new WLFTransport<Thermo>(*this)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::autoPtr<Foam::WLFTransport<Thermo>>
|
||||||
|
Foam::WLFTransport<Thermo>::New
|
||||||
|
(
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
{
|
||||||
|
return autoPtr<WLFTransport<Thermo>>
|
||||||
|
(
|
||||||
|
new WLFTransport<Thermo>(dict)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::scalar Foam::WLFTransport<Thermo>::mu
|
||||||
|
(
|
||||||
|
const scalar p,
|
||||||
|
const scalar T
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
return mu0_*exp(-C1_*(T - Tr_)/(C2_ + T - Tr_));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::scalar Foam::WLFTransport<Thermo>::kappa
|
||||||
|
(
|
||||||
|
const scalar p, const scalar T
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
return this->Cp(p, T)*mu(p, T)*rPr_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::scalar Foam::WLFTransport<Thermo>::alphah
|
||||||
|
(
|
||||||
|
const scalar p,
|
||||||
|
const scalar T
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
|
||||||
|
return mu(p, T)*rPr_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline void Foam::WLFTransport<Thermo>::operator=
|
||||||
|
(
|
||||||
|
const WLFTransport<Thermo>& wlft
|
||||||
|
)
|
||||||
|
{
|
||||||
|
Thermo::operator=(wlft);
|
||||||
|
|
||||||
|
mu0_ = wlft.mu0_;
|
||||||
|
Tr_ = wlft.Tr_;
|
||||||
|
C1_ = wlft.C1_;
|
||||||
|
C2_ = wlft.C2_;
|
||||||
|
rPr_ = wlft.rPr_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline void Foam::WLFTransport<Thermo>::operator+=
|
||||||
|
(
|
||||||
|
const WLFTransport<Thermo>& wlft
|
||||||
|
)
|
||||||
|
{
|
||||||
|
scalar Y1 = this->Y();
|
||||||
|
|
||||||
|
Thermo::operator+=(wlft);
|
||||||
|
|
||||||
|
if (mag(this->Y()) > small)
|
||||||
|
{
|
||||||
|
Y1 /= this->Y();
|
||||||
|
scalar Y2 = wlft.Y()/this->Y();
|
||||||
|
|
||||||
|
mu0_ = Y1*mu0_ + Y2*wlft.mu0_;
|
||||||
|
Tr_ = Y1*Tr_ + Y2*wlft.Tr_;
|
||||||
|
C1_ = Y1*C1_ + Y2*wlft.C1_;
|
||||||
|
C2_ = Y1*C2_ + Y2*wlft.C2_;
|
||||||
|
rPr_ = 1.0/(Y1/rPr_ + Y2/wlft.rPr_);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline void Foam::WLFTransport<Thermo>::operator*=
|
||||||
|
(
|
||||||
|
const scalar s
|
||||||
|
)
|
||||||
|
{
|
||||||
|
Thermo::operator*=(s);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::WLFTransport<Thermo> Foam::operator+
|
||||||
|
(
|
||||||
|
const WLFTransport<Thermo>& wlft1,
|
||||||
|
const WLFTransport<Thermo>& wlft2
|
||||||
|
)
|
||||||
|
{
|
||||||
|
Thermo t
|
||||||
|
(
|
||||||
|
static_cast<const Thermo&>(wlft1) + static_cast<const Thermo&>(wlft2)
|
||||||
|
);
|
||||||
|
|
||||||
|
if (mag(t.Y()) < small)
|
||||||
|
{
|
||||||
|
return WLFTransport<Thermo>
|
||||||
|
(
|
||||||
|
t,
|
||||||
|
0,
|
||||||
|
wlft1.mu0_,
|
||||||
|
wlft1.Tr_,
|
||||||
|
wlft1.C1_,
|
||||||
|
wlft1.C2_,
|
||||||
|
wlft1.rPr_
|
||||||
|
);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
scalar Y1 = wlft1.Y()/t.Y();
|
||||||
|
scalar Y2 = wlft2.Y()/t.Y();
|
||||||
|
|
||||||
|
return WLFTransport<Thermo>
|
||||||
|
(
|
||||||
|
t,
|
||||||
|
Y1*wlft1.mu0_ + Y2*wlft2.mu0_,
|
||||||
|
Y1*wlft1.Tr_ + Y2*wlft2.Tr_,
|
||||||
|
Y1*wlft1.C1_ + Y2*wlft2.C1_,
|
||||||
|
Y1*wlft1.C2_ + Y2*wlft2.C2_,
|
||||||
|
1.0/(Y1/wlft1.rPr_ + Y2/wlft2.rPr_)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class Thermo>
|
||||||
|
inline Foam::WLFTransport<Thermo> Foam::operator*
|
||||||
|
(
|
||||||
|
const scalar s,
|
||||||
|
const WLFTransport<Thermo>& wlft
|
||||||
|
)
|
||||||
|
{
|
||||||
|
return WLFTransport<Thermo>
|
||||||
|
(
|
||||||
|
s*static_cast<const Thermo&>(wlft),
|
||||||
|
wlft.mu0_,
|
||||||
|
wlft.Tr_,
|
||||||
|
wlft.C1_,
|
||||||
|
wlft.C2_,
|
||||||
|
1.0/wlft.rPr_
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
Reference in New Issue
Block a user