Skip to content
Snippets Groups Projects
Commit a77070cd authored by Henry Weller's avatar Henry Weller
Browse files

kOmegaSSTSato: Relocated to src/TurbulenceModels/phaseCompressible/RAS

parent 2127051a
Branches
Tags
No related merge requests found
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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 "kOmegaSSTSato.H"
#include "addToRunTimeSelectionTable.H"
#include "twoPhaseSystem.H"
#include "dragModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace RASModels
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
kOmegaSSTSato<BasicTurbulenceModel>::kOmegaSSTSato
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName,
const word& type
)
:
kOmegaSST<BasicTurbulenceModel>
(
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport,
propertiesName,
type
),
gasTurbulencePtr_(NULL),
Cmub_
(
dimensioned<scalar>::lookupOrAddToDict
(
"Cmub",
this->coeffDict_,
0.6
)
)
{
if (type == typeName)
{
correctNut();
this->printCoeffs(type);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
bool kOmegaSSTSato<BasicTurbulenceModel>::read()
{
if (kOmegaSST<BasicTurbulenceModel>::read())
{
Cmub_.readIfPresent(this->coeffDict());
return true;
}
else
{
return false;
}
}
template<class BasicTurbulenceModel>
const PhaseCompressibleTurbulenceModel
<
typename BasicTurbulenceModel::transportModel
>&
kOmegaSSTSato<BasicTurbulenceModel>::gasTurbulence() const
{
if (!gasTurbulencePtr_)
{
const volVectorField& U = this->U_;
const transportModel& liquid = this->transport();
const twoPhaseSystem& fluid = liquid.fluid();
const transportModel& gas = fluid.otherPhase(liquid);
gasTurbulencePtr_ =
&U.db()
.lookupObject<PhaseCompressibleTurbulenceModel<transportModel> >
(
IOobject::groupName
(
turbulenceModel::propertiesName,
gas.name()
)
);
}
return *gasTurbulencePtr_;
}
template<class BasicTurbulenceModel>
void kOmegaSSTSato<BasicTurbulenceModel>::correctNut()
{
const PhaseCompressibleTurbulenceModel<transportModel>& gasTurbulence =
this->gasTurbulence();
volScalarField yPlus
(
pow(this->betaStar_, 0.25)*this->y_*sqrt(this->k_)/this->nu()
);
this->nut_ =
this->a1_*this->k_
/max
(
this->a1_*this->omega_,
this->F23()*sqrt(2.0)*mag(symm(fvc::grad(this->U_)))
)
+ sqr(1 - exp(-yPlus/16.0))
*Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
*(mag(this->U_ - gasTurbulence.U()));
this->nut_.correctBoundaryConditions();
}
template<class BasicTurbulenceModel>
void kOmegaSSTSato<BasicTurbulenceModel>::correct()
{
kOmegaSST<BasicTurbulenceModel>::correct();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace Foam
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014-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::RASModels::kOmegaSSTSato
Group
grpRASTurbulence
Description
Implementation of the k-omega-SST turbulence model for dispersed bubbly
flows with Sato (1981) bubble induced turbulent viscosity model.
Bubble induced turbulent viscosity model described in:
\verbatim
Sato, Y., Sadatomi, M.
"Momentum and heat transfer in two-phase bubble flow - I, Theory"
International Journal of Multiphase FLow 7, pp. 167-177, 1981.
\endverbatim
Turbulence model described in:
\verbatim
Menter, F., Esch, T.
"Elements of Industrial Heat Transfer Prediction"
16th Brazilian Congress of Mechanical Engineering (COBEM),
Nov. 2001
\endverbatim
with the addition of the optional F3 term for rough walls from
\verbatim
Hellsten, A.
"Some Improvements in Menter’s k-omega-SST turbulence model"
29th AIAA Fluid Dynamics Conference,
AIAA-98-2554,
June 1998.
\endverbatim
Note that this implementation is written in terms of alpha diffusion
coefficients rather than the more traditional sigma (alpha = 1/sigma) so
that the blending can be applied to all coefficuients in a consistent
manner. The paper suggests that sigma is blended but this would not be
consistent with the blending of the k-epsilon and k-omega models.
Also note that the error in the last term of equation (2) relating to
sigma has been corrected.
Wall-functions are applied in this implementation by using equations (14)
to specify the near-wall omega as appropriate.
The blending functions (15) and (16) are not currently used because of the
uncertainty in their origin, range of applicability and that is y+ becomes
sufficiently small blending u_tau in this manner clearly becomes nonsense.
The default model coefficients correspond to the following:
\verbatim
kOmegaSSTCoeffs
{
alphaK1 0.85034;
alphaK2 1.0;
alphaOmega1 0.5;
alphaOmega2 0.85616;
Prt 1.0; // only for compressible
beta1 0.075;
beta2 0.0828;
betaStar 0.09;
gamma1 0.5532;
gamma2 0.4403;
a1 0.31;
b1 1.0;
c1 10.0;
F3 no;
Cmub 0.6;
}
\endverbatim
SourceFiles
kOmegaSST.C
SourceFiles
kOmegaSSTSato.C
\*---------------------------------------------------------------------------*/
#ifndef kOmegaSSTSato_H
#define kOmegaSSTSato_H
#include "kOmegaSST.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class kOmegaSSTSato Declaration
\*---------------------------------------------------------------------------*/
template<class BasicTurbulenceModel>
class kOmegaSSTSato
:
public kOmegaSST<BasicTurbulenceModel>
{
// Private data
mutable const PhaseCompressibleTurbulenceModel
<
typename BasicTurbulenceModel::transportModel
> *gasTurbulencePtr_;
// Private Member Functions
//- Return the turbulence model for the gas phase
const PhaseCompressibleTurbulenceModel
<
typename BasicTurbulenceModel::transportModel
>&
gasTurbulence() const;
// Disallow default bitwise copy construct and assignment
kOmegaSSTSato(const kOmegaSSTSato&);
kOmegaSSTSato& operator=(const kOmegaSSTSato&);
protected:
// Protected data
// Model coefficients
dimensionedScalar Cmub_;
// Protected Member Functions
virtual void correctNut();
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("kOmegaSSTSato");
// Constructors
//- Construct from components
kOmegaSSTSato
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName = turbulenceModel::propertiesName,
const word& type = typeName
);
//- Destructor
virtual ~kOmegaSSTSato()
{}
// Member Functions
//- Read model coefficients if they have changed
virtual bool read();
//- Solve the turbulence equations and correct the turbulence viscosity
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "kOmegaSSTSato.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment