mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
equationOfState/Boussinesq: New equation of state for the Boussinesq approximation for buoyant flows
Description
Incompressible gas equation of state using the Boussinesq approximation for
the density as a function of temperature only:
\verbatim
rho = rho0*(1 - beta*(T - T0))
\endverbatim
To be used with the buoyantPimpleFoam and buoyantSimpleFoam solvers as
an alternative to using buoyantBoussinesqPimpleFoam or
buoyantBoussinesqSimpleFoam, providing consistency with all other
solvers and utilities using the thermodynamics package in OpenFOAM.
This commit is contained in:
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -29,6 +29,7 @@ License
|
||||
#include "specie.H"
|
||||
#include "perfectGas.H"
|
||||
#include "incompressiblePerfectGas.H"
|
||||
#include "Boussinesq.H"
|
||||
#include "rhoConst.H"
|
||||
#include "perfectFluid.H"
|
||||
#include "PengRobinsonGas.H"
|
||||
@ -177,6 +178,42 @@ makeThermo
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
constTransport,
|
||||
sensibleEnthalpy,
|
||||
hConstThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
sutherlandTransport,
|
||||
sensibleEnthalpy,
|
||||
hConstThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
sutherlandTransport,
|
||||
sensibleEnthalpy,
|
||||
janafThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
@ -335,6 +372,42 @@ makeThermo
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
constTransport,
|
||||
sensibleInternalEnergy,
|
||||
hConstThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
sutherlandTransport,
|
||||
sensibleInternalEnergy,
|
||||
hConstThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
makeThermo
|
||||
(
|
||||
rhoThermo,
|
||||
heRhoThermo,
|
||||
pureMixture,
|
||||
sutherlandTransport,
|
||||
sensibleInternalEnergy,
|
||||
janafThermo,
|
||||
Boussinesq,
|
||||
specie
|
||||
);
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
|
||||
@ -0,0 +1,98 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2015 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 "Boussinesq.H"
|
||||
#include "IOstreams.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
Foam::Boussinesq<Specie>::Boussinesq(Istream& is)
|
||||
:
|
||||
Specie(is),
|
||||
rho0_(readScalar(is)),
|
||||
T0_(readScalar(is)),
|
||||
beta_(readScalar(is))
|
||||
{
|
||||
is.check
|
||||
(
|
||||
"Boussinesq<Specie>::"
|
||||
"Boussinesq(Istream& is)"
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
Foam::Boussinesq<Specie>::Boussinesq
|
||||
(
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
Specie(dict),
|
||||
rho0_(readScalar(dict.subDict("equationOfState").lookup("rho0"))),
|
||||
T0_(readScalar(dict.subDict("equationOfState").lookup("T0"))),
|
||||
beta_(readScalar(dict.subDict("equationOfState").lookup("beta")))
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
void Foam::Boussinesq<Specie>::write(Ostream& os) const
|
||||
{
|
||||
Specie::write(os);
|
||||
dictionary dict("equationOfState");
|
||||
dict.add("rho0", rho0_);
|
||||
dict.add("T0", T0_);
|
||||
dict.add("beta", beta_);
|
||||
|
||||
os << indent << dict.dictName() << dict;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
Foam::Ostream& Foam::operator<<
|
||||
(
|
||||
Ostream& os,
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
{
|
||||
os << static_cast<const Specie&>(b)
|
||||
<< token::SPACE << b.rho0_
|
||||
<< token::SPACE << b.T0_
|
||||
<< token::SPACE << b.beta_;
|
||||
|
||||
os.check
|
||||
(
|
||||
"Ostream& operator<<"
|
||||
"(Ostream& os, const Boussinesq<Specie>& st)"
|
||||
);
|
||||
return os;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,260 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2015 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::Boussinesq
|
||||
|
||||
Description
|
||||
Incompressible gas equation of state using the Boussinesq approximation for
|
||||
the density as a function of temperature only:
|
||||
|
||||
\verbatim
|
||||
rho = rho0*(1 - beta*(T - T0))
|
||||
\endverbatim
|
||||
|
||||
SourceFiles
|
||||
BoussinesqI.H
|
||||
Boussinesq.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef Boussinesq_H
|
||||
#define Boussinesq_H
|
||||
|
||||
#include "autoPtr.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
// Forward declaration of friend functions and operators
|
||||
|
||||
template<class Specie> class Boussinesq;
|
||||
|
||||
template<class Specie>
|
||||
inline Boussinesq<Specie> operator+
|
||||
(
|
||||
const Boussinesq<Specie>&,
|
||||
const Boussinesq<Specie>&
|
||||
);
|
||||
|
||||
template<class Specie>
|
||||
inline Boussinesq<Specie> operator-
|
||||
(
|
||||
const Boussinesq<Specie>&,
|
||||
const Boussinesq<Specie>&
|
||||
);
|
||||
|
||||
template<class Specie>
|
||||
inline Boussinesq<Specie> operator*
|
||||
(
|
||||
const scalar,
|
||||
const Boussinesq<Specie>&
|
||||
);
|
||||
|
||||
template<class Specie>
|
||||
inline Boussinesq<Specie> operator==
|
||||
(
|
||||
const Boussinesq<Specie>&,
|
||||
const Boussinesq<Specie>&
|
||||
);
|
||||
|
||||
template<class Specie>
|
||||
Ostream& operator<<
|
||||
(
|
||||
Ostream&,
|
||||
const Boussinesq<Specie>&
|
||||
);
|
||||
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class Boussinesq Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
template<class Specie>
|
||||
class Boussinesq
|
||||
:
|
||||
public Specie
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- Reference density
|
||||
scalar rho0_;
|
||||
|
||||
//- Reference temperature
|
||||
scalar T0_;
|
||||
|
||||
//- Thermal expansion coefficient
|
||||
scalar beta_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
inline Boussinesq
|
||||
(
|
||||
const Specie& sp,
|
||||
const scalar rho0,
|
||||
const scalar T0,
|
||||
const scalar beta
|
||||
);
|
||||
|
||||
//- Construct from Boussinesq
|
||||
inline Boussinesq(const Boussinesq& sp);
|
||||
|
||||
//- Construct from Istream
|
||||
Boussinesq(Istream&);
|
||||
|
||||
//- Construct from dictionary
|
||||
Boussinesq(const dictionary& dict);
|
||||
|
||||
//- Construct as named copy
|
||||
inline Boussinesq
|
||||
(
|
||||
const word& name,
|
||||
const Boussinesq&
|
||||
);
|
||||
|
||||
//- Construct and return a clone
|
||||
inline autoPtr<Boussinesq> clone() const;
|
||||
|
||||
// Selector from Istream
|
||||
inline static autoPtr<Boussinesq> New(Istream& is);
|
||||
|
||||
// Selector from dictionary
|
||||
inline static autoPtr<Boussinesq> New
|
||||
(
|
||||
const dictionary& dict
|
||||
);
|
||||
|
||||
|
||||
// Member functions
|
||||
|
||||
//- Return the instantiated type name
|
||||
static word typeName()
|
||||
{
|
||||
return
|
||||
"Boussinesq<"
|
||||
+ word(Specie::typeName_()) + '>';
|
||||
}
|
||||
|
||||
|
||||
// Fundamental properties
|
||||
|
||||
//- Is the equation of state is incompressible i.e. rho != f(p)
|
||||
static const bool incompressible = true;
|
||||
|
||||
//- Is the equation of state is isochoric i.e. rho = const
|
||||
static const bool isochoric = false;
|
||||
|
||||
//- Return density [kg/m^3]
|
||||
inline scalar rho(scalar p, scalar T) const;
|
||||
|
||||
//- Return entropy [J/(kmol K)]
|
||||
inline scalar s(const scalar p, const scalar T) const;
|
||||
|
||||
//- Return compressibility rho/p [s^2/m^2]
|
||||
inline scalar psi(scalar p, scalar T) const;
|
||||
|
||||
//- Return compression factor []
|
||||
inline scalar Z(scalar p, scalar T) const;
|
||||
|
||||
//- Return (cp - cv) [J/(kmol K]
|
||||
inline scalar cpMcv(scalar p, scalar T) const;
|
||||
|
||||
|
||||
// IO
|
||||
|
||||
//- Write to Ostream
|
||||
void write(Ostream& os) const;
|
||||
|
||||
|
||||
// Member operators
|
||||
|
||||
inline Boussinesq& operator=
|
||||
(
|
||||
const Boussinesq&
|
||||
);
|
||||
inline void operator+=(const Boussinesq&);
|
||||
inline void operator-=(const Boussinesq&);
|
||||
|
||||
inline void operator*=(const scalar);
|
||||
|
||||
|
||||
// Friend operators
|
||||
|
||||
friend Boussinesq operator+ <Specie>
|
||||
(
|
||||
const Boussinesq&,
|
||||
const Boussinesq&
|
||||
);
|
||||
|
||||
friend Boussinesq operator- <Specie>
|
||||
(
|
||||
const Boussinesq&,
|
||||
const Boussinesq&
|
||||
);
|
||||
|
||||
friend Boussinesq operator* <Specie>
|
||||
(
|
||||
const scalar s,
|
||||
const Boussinesq&
|
||||
);
|
||||
|
||||
friend Boussinesq operator== <Specie>
|
||||
(
|
||||
const Boussinesq&,
|
||||
const Boussinesq&
|
||||
);
|
||||
|
||||
|
||||
// Ostream Operator
|
||||
|
||||
friend Ostream& operator<< <Specie>
|
||||
(
|
||||
Ostream&,
|
||||
const Boussinesq&
|
||||
);
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#include "BoussinesqI.H"
|
||||
|
||||
#ifdef NoRepository
|
||||
# include "Boussinesq.C"
|
||||
#endif
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,292 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2015 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 "Boussinesq.H"
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie>::Boussinesq
|
||||
(
|
||||
const Specie& sp, const scalar rho0, const scalar T0, const scalar beta
|
||||
)
|
||||
:
|
||||
Specie(sp),
|
||||
rho0_(rho0),
|
||||
T0_(T0),
|
||||
beta_(beta)
|
||||
{}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie>::Boussinesq
|
||||
(
|
||||
const Boussinesq& b
|
||||
)
|
||||
:
|
||||
Specie(b),
|
||||
rho0_(b.rho0_),
|
||||
T0_(b.T0_),
|
||||
beta_(b.beta_)
|
||||
{}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie>::Boussinesq
|
||||
(
|
||||
const word& name,
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
:
|
||||
Specie(name, b),
|
||||
rho0_(b.rho0_),
|
||||
T0_(b.T0_),
|
||||
beta_(b.beta_)
|
||||
{}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::autoPtr<Foam::Boussinesq<Specie> >
|
||||
Foam::Boussinesq<Specie>::clone() const
|
||||
{
|
||||
return autoPtr<Boussinesq<Specie> >
|
||||
(
|
||||
new Boussinesq<Specie>(*this)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::autoPtr<Foam::Boussinesq<Specie> >
|
||||
Foam::Boussinesq<Specie>::New
|
||||
(
|
||||
Istream& is
|
||||
)
|
||||
{
|
||||
return autoPtr<Boussinesq<Specie> >
|
||||
(
|
||||
new Boussinesq<Specie>(is)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::autoPtr<Foam::Boussinesq<Specie> >
|
||||
Foam::Boussinesq<Specie>::New
|
||||
(
|
||||
const dictionary& dict
|
||||
)
|
||||
{
|
||||
return autoPtr<Boussinesq<Specie> >
|
||||
(
|
||||
new Boussinesq<Specie>(dict)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::scalar Foam::Boussinesq<Specie>::rho
|
||||
(
|
||||
scalar p,
|
||||
scalar T
|
||||
) const
|
||||
{
|
||||
return rho0_*(1.0 - beta_*(T - T0_));
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::scalar Foam::Boussinesq<Specie>::s
|
||||
(
|
||||
scalar p,
|
||||
scalar T
|
||||
) const
|
||||
{
|
||||
return 0;
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::scalar Foam::Boussinesq<Specie>::psi
|
||||
(
|
||||
scalar p,
|
||||
scalar T
|
||||
) const
|
||||
{
|
||||
return 0;
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::scalar Foam::Boussinesq<Specie>::Z
|
||||
(
|
||||
scalar p,
|
||||
scalar T
|
||||
) const
|
||||
{
|
||||
return 0;
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::scalar Foam::Boussinesq<Specie>::cpMcv
|
||||
(
|
||||
scalar p,
|
||||
scalar T
|
||||
) const
|
||||
{
|
||||
return RR;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie>&
|
||||
Foam::Boussinesq<Specie>::operator=
|
||||
(
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
{
|
||||
Specie::operator=(b);
|
||||
|
||||
rho0_ = b.rho0_;
|
||||
T0_ = b.T0_;
|
||||
beta_ = b.beta_;
|
||||
|
||||
return *this;
|
||||
}
|
||||
|
||||
template<class Specie>
|
||||
inline void Foam::Boussinesq<Specie>::operator+=
|
||||
(
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
{
|
||||
scalar molr1 = this->nMoles();
|
||||
Specie::operator+=(b);
|
||||
molr1 /= this->nMoles();
|
||||
scalar molr2 = b.nMoles()/this->nMoles();
|
||||
|
||||
rho0_ = molr1*rho0_ + molr2*b.rho0_;
|
||||
T0_ = molr1*T0_ + molr2*b.T0_;
|
||||
beta_ = molr1*beta_ + molr2*b.beta_;
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline void Foam::Boussinesq<Specie>::operator-=
|
||||
(
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
{
|
||||
Specie::operator-=(b);
|
||||
rho0_ = b.rho0_;
|
||||
T0_ = b.T0_;
|
||||
beta_ = b.beta_;
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline void Foam::Boussinesq<Specie>::operator*=(const scalar s)
|
||||
{
|
||||
Specie::operator*=(s);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie> Foam::operator+
|
||||
(
|
||||
const Boussinesq<Specie>& b1,
|
||||
const Boussinesq<Specie>& b2
|
||||
)
|
||||
{
|
||||
scalar nMoles = b1.nMoles() + b2.nMoles();
|
||||
scalar molr1 = b1.nMoles()/nMoles;
|
||||
scalar molr2 = b2.nMoles()/nMoles;
|
||||
|
||||
return Boussinesq<Specie>
|
||||
(
|
||||
static_cast<const Specie&>(b1)
|
||||
+ static_cast<const Specie&>(b2),
|
||||
molr1*b1.rho0_ + molr2*b2.rho0_,
|
||||
molr1*b1.T0_ + molr2*b2.T0_,
|
||||
molr1*b1.beta_ + molr2*b2.beta_
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie> Foam::operator-
|
||||
(
|
||||
const Boussinesq<Specie>& b1,
|
||||
const Boussinesq<Specie>& b2
|
||||
)
|
||||
{
|
||||
return Boussinesq<Specie>
|
||||
(
|
||||
static_cast<const Specie&>(b1)
|
||||
- static_cast<const Specie&>(b2),
|
||||
b1.rho0_ - b2.rho0_,
|
||||
b1.T0_ - b2.T0_,
|
||||
b1.beta_ - b2.beta_
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie> Foam::operator*
|
||||
(
|
||||
const scalar s,
|
||||
const Boussinesq<Specie>& b
|
||||
)
|
||||
{
|
||||
return Boussinesq<Specie>
|
||||
(
|
||||
s*static_cast<const Specie&>(b),
|
||||
b.rho0_,
|
||||
b.T0_,
|
||||
b.beta_
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class Specie>
|
||||
inline Foam::Boussinesq<Specie> Foam::operator==
|
||||
(
|
||||
const Boussinesq<Specie>& pg1,
|
||||
const Boussinesq<Specie>& pg2
|
||||
)
|
||||
{
|
||||
return pg2 - pg1;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
Reference in New Issue
Block a user