mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
turbulenceModels/LES: Added WALE model
Changed the tutorials/incompressible/pimpleFoam/channel395 to demonstrate the WALE model
This commit is contained in:
@ -72,6 +72,9 @@ makeRASModel(kOmegaSST);
|
|||||||
#include "Smagorinsky.H"
|
#include "Smagorinsky.H"
|
||||||
makeLESModel(Smagorinsky);
|
makeLESModel(Smagorinsky);
|
||||||
|
|
||||||
|
#include "WALE.H"
|
||||||
|
makeLESModel(WALE);
|
||||||
|
|
||||||
#include "kEqn.H"
|
#include "kEqn.H"
|
||||||
makeLESModel(kEqn);
|
makeLESModel(kEqn);
|
||||||
|
|
||||||
|
|||||||
@ -67,6 +67,9 @@ makeRASModel(kOmegaSST);
|
|||||||
#include "Smagorinsky.H"
|
#include "Smagorinsky.H"
|
||||||
makeLESModel(Smagorinsky);
|
makeLESModel(Smagorinsky);
|
||||||
|
|
||||||
|
#include "WALE.H"
|
||||||
|
makeLESModel(WALE);
|
||||||
|
|
||||||
#include "kEqn.H"
|
#include "kEqn.H"
|
||||||
makeLESModel(kEqn);
|
makeLESModel(kEqn);
|
||||||
|
|
||||||
|
|||||||
205
src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.C
Normal file
205
src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.C
Normal file
@ -0,0 +1,205 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / 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 "WALE.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace LESModels
|
||||||
|
{
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
tmp<volSymmTensorField> WALE<BasicTurbulenceModel>::Sd
|
||||||
|
(
|
||||||
|
const volTensorField& gradU
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
return dev(symm(gradU & gradU));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
tmp<volScalarField> WALE<BasicTurbulenceModel>::k
|
||||||
|
(
|
||||||
|
const volTensorField& gradU
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
volScalarField magSqrSd(magSqr(Sd(gradU)));
|
||||||
|
|
||||||
|
return tmp<volScalarField>
|
||||||
|
(
|
||||||
|
new volScalarField
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
IOobject::groupName("k", this->U_.group()),
|
||||||
|
this->runTime_.timeName(),
|
||||||
|
this->mesh_
|
||||||
|
),
|
||||||
|
sqr(sqr(Cw_)*this->delta()/Ck_)*
|
||||||
|
(
|
||||||
|
pow3(magSqrSd)
|
||||||
|
/(
|
||||||
|
sqr
|
||||||
|
(
|
||||||
|
pow(magSqr(symm(gradU)), 5.0/2.0)
|
||||||
|
+ pow(magSqrSd, 5.0/4.0)
|
||||||
|
)
|
||||||
|
+ dimensionedScalar
|
||||||
|
(
|
||||||
|
"SMALL",
|
||||||
|
dimensionSet(0, 0, -10, 0, 0),
|
||||||
|
SMALL
|
||||||
|
)
|
||||||
|
)
|
||||||
|
)
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
void WALE<BasicTurbulenceModel>::correctNut()
|
||||||
|
{
|
||||||
|
this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
|
||||||
|
this->nut_.correctBoundaryConditions();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
WALE<BasicTurbulenceModel>::WALE
|
||||||
|
(
|
||||||
|
const alphaField& alpha,
|
||||||
|
const rhoField& rho,
|
||||||
|
const volVectorField& U,
|
||||||
|
const surfaceScalarField& alphaRhoPhi,
|
||||||
|
const surfaceScalarField& phi,
|
||||||
|
const transportModel& transport,
|
||||||
|
const word& propertiesName,
|
||||||
|
const word& type
|
||||||
|
)
|
||||||
|
:
|
||||||
|
LESeddyViscosity<BasicTurbulenceModel>
|
||||||
|
(
|
||||||
|
type,
|
||||||
|
alpha,
|
||||||
|
rho,
|
||||||
|
U,
|
||||||
|
alphaRhoPhi,
|
||||||
|
phi,
|
||||||
|
transport,
|
||||||
|
propertiesName
|
||||||
|
),
|
||||||
|
|
||||||
|
Ck_
|
||||||
|
(
|
||||||
|
dimensioned<scalar>::lookupOrAddToDict
|
||||||
|
(
|
||||||
|
"Ck",
|
||||||
|
this->coeffDict_,
|
||||||
|
0.094
|
||||||
|
)
|
||||||
|
),
|
||||||
|
|
||||||
|
Cw_
|
||||||
|
(
|
||||||
|
dimensioned<scalar>::lookupOrAddToDict
|
||||||
|
(
|
||||||
|
"Cw",
|
||||||
|
this->coeffDict_,
|
||||||
|
0.325
|
||||||
|
)
|
||||||
|
)
|
||||||
|
{
|
||||||
|
if (type == typeName)
|
||||||
|
{
|
||||||
|
correctNut();
|
||||||
|
this->printCoeffs(type);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
bool WALE<BasicTurbulenceModel>::read()
|
||||||
|
{
|
||||||
|
if (LESeddyViscosity<BasicTurbulenceModel>::read())
|
||||||
|
{
|
||||||
|
Ck_.readIfPresent(this->coeffDict());
|
||||||
|
Cw_.readIfPresent(this->coeffDict());
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
tmp<volScalarField> WALE<BasicTurbulenceModel>::epsilon() const
|
||||||
|
{
|
||||||
|
volScalarField k(this->k(fvc::grad(this->U_)));
|
||||||
|
|
||||||
|
return tmp<volScalarField>
|
||||||
|
(
|
||||||
|
new volScalarField
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
IOobject::groupName("epsilon", this->U_.group()),
|
||||||
|
this->runTime_.timeName(),
|
||||||
|
this->mesh_,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE
|
||||||
|
),
|
||||||
|
this->Ce_*k*sqrt(k)/this->delta()
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
void WALE<BasicTurbulenceModel>::correct()
|
||||||
|
{
|
||||||
|
LESeddyViscosity<BasicTurbulenceModel>::correct();
|
||||||
|
correctNut();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace LESModels
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
173
src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.H
Normal file
173
src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.H
Normal file
@ -0,0 +1,173 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / 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::LESModels::WALE
|
||||||
|
|
||||||
|
Group
|
||||||
|
grpLESTurbulence
|
||||||
|
|
||||||
|
Description
|
||||||
|
The Wall-adapting local eddy-viscosity (WALE) SGS model.
|
||||||
|
|
||||||
|
Reference:
|
||||||
|
\verbatim
|
||||||
|
Nicoud, F., & Ducros, F. (1999).
|
||||||
|
Subgrid-scale stress modelling based on the square of the velocity
|
||||||
|
gradient tensor.
|
||||||
|
Flow, Turbulence and Combustion, 62(3), 183-200.
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
The default model coefficients correspond to the following:
|
||||||
|
\verbatim
|
||||||
|
WALECoeffs
|
||||||
|
{
|
||||||
|
Ck 0.094;
|
||||||
|
Ce 1.048;
|
||||||
|
Cw 0.325;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
WALE.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef WALE_H
|
||||||
|
#define WALE_H
|
||||||
|
|
||||||
|
#include "LESModel.H"
|
||||||
|
#include "LESeddyViscosity.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace LESModels
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class WALE Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
template<class BasicTurbulenceModel>
|
||||||
|
class WALE
|
||||||
|
:
|
||||||
|
public LESeddyViscosity<BasicTurbulenceModel>
|
||||||
|
{
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
// Disallow default bitwise copy construct and assignment
|
||||||
|
WALE(const WALE&);
|
||||||
|
WALE& operator=(const WALE&);
|
||||||
|
|
||||||
|
|
||||||
|
protected:
|
||||||
|
|
||||||
|
// Protected data
|
||||||
|
|
||||||
|
dimensionedScalar Ck_;
|
||||||
|
dimensionedScalar Cw_;
|
||||||
|
|
||||||
|
|
||||||
|
// Protected Member Functions
|
||||||
|
|
||||||
|
//- Return the deviatoric symmetric part of the square of the given
|
||||||
|
// velocity gradient field
|
||||||
|
tmp<volSymmTensorField> Sd(const volTensorField& gradU) const;
|
||||||
|
|
||||||
|
//- Return SGS kinetic energy
|
||||||
|
// calculated from the given velocity gradient
|
||||||
|
tmp<volScalarField> k(const volTensorField& gradU) const;
|
||||||
|
|
||||||
|
//- Update the SGS eddy viscosity
|
||||||
|
virtual void correctNut();
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
typedef typename BasicTurbulenceModel::alphaField alphaField;
|
||||||
|
typedef typename BasicTurbulenceModel::rhoField rhoField;
|
||||||
|
typedef typename BasicTurbulenceModel::transportModel transportModel;
|
||||||
|
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("WALE");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
WALE
|
||||||
|
(
|
||||||
|
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 ~WALE()
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Read model coefficients if they have changed
|
||||||
|
virtual bool read();
|
||||||
|
|
||||||
|
//- Return SGS kinetic energy
|
||||||
|
virtual tmp<volScalarField> k() const
|
||||||
|
{
|
||||||
|
return k(fvc::grad(this->U_));
|
||||||
|
}
|
||||||
|
|
||||||
|
//- Return sub-grid disipation rate
|
||||||
|
virtual tmp<volScalarField> epsilon() const;
|
||||||
|
|
||||||
|
//- Correct Eddy-Viscosity and related properties
|
||||||
|
virtual void correct();
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace LESModels
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#ifdef NoRepository
|
||||||
|
# include "WALE.C"
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -19,13 +19,13 @@ simulationType LES;
|
|||||||
|
|
||||||
LES
|
LES
|
||||||
{
|
{
|
||||||
LESModel SpalartAllmarasDDES; //kEqn;
|
LESModel WALE;
|
||||||
|
|
||||||
turbulence on;
|
turbulence on;
|
||||||
|
|
||||||
printCoeffs on;
|
printCoeffs on;
|
||||||
|
|
||||||
delta vanDriest;
|
delta cubeRootVol;
|
||||||
|
|
||||||
cubeRootVolCoeffs
|
cubeRootVolCoeffs
|
||||||
{
|
{
|
||||||
|
|||||||
Reference in New Issue
Block a user