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

turbulenceModels/LES/dynamicLagrangian: Added the Lagrangian averaged form of...

turbulenceModels/LES/dynamicLagrangian: Added the Lagrangian averaged form of the dynamic Smagorinsky eddy-viscosity SGS model for LES
parent 8628ef2f
No related branches found
No related tags found
No related merge requests found
...@@ -34,6 +34,8 @@ License ...@@ -34,6 +34,8 @@ License
#include "RASModel.H" #include "RASModel.H"
#include "LESModel.H" #include "LESModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makeBaseTurbulenceModel makeBaseTurbulenceModel
( (
geometricOneField, geometricOneField,
...@@ -51,6 +53,11 @@ makeBaseTurbulenceModel ...@@ -51,6 +53,11 @@ makeBaseTurbulenceModel
makeTemplatedTurbulenceModel \ makeTemplatedTurbulenceModel \
(fluidThermothermalDiffusivity, LES, Type) (fluidThermothermalDiffusivity, LES, Type)
// -------------------------------------------------------------------------- //
// RAS models
// -------------------------------------------------------------------------- //
#include "SpalartAllmaras.H" #include "SpalartAllmaras.H"
makeRASModel(SpalartAllmaras); makeRASModel(SpalartAllmaras);
...@@ -85,12 +92,19 @@ makeRASModel(LRR); ...@@ -85,12 +92,19 @@ makeRASModel(LRR);
makeRASModel(SSG); makeRASModel(SSG);
// -------------------------------------------------------------------------- //
// LES models
// -------------------------------------------------------------------------- //
#include "Smagorinsky.H" #include "Smagorinsky.H"
makeLESModel(Smagorinsky); makeLESModel(Smagorinsky);
#include "WALE.H" #include "WALE.H"
makeLESModel(WALE); makeLESModel(WALE);
#include "dynamicLagrangian.H"
makeLESModel(dynamicLagrangian);
#include "kEqn.H" #include "kEqn.H"
makeLESModel(kEqn); makeLESModel(kEqn);
......
...@@ -32,6 +32,8 @@ License ...@@ -32,6 +32,8 @@ License
#include "RASModel.H" #include "RASModel.H"
#include "LESModel.H" #include "LESModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makeBaseTurbulenceModel makeBaseTurbulenceModel
( (
geometricOneField, geometricOneField,
...@@ -49,6 +51,11 @@ makeBaseTurbulenceModel ...@@ -49,6 +51,11 @@ makeBaseTurbulenceModel
makeTemplatedTurbulenceModel \ makeTemplatedTurbulenceModel \
(transportModelIncompressibleTurbulenceModel, LES, Type) (transportModelIncompressibleTurbulenceModel, LES, Type)
// -------------------------------------------------------------------------- //
// RAS models
// -------------------------------------------------------------------------- //
#include "SpalartAllmaras.H" #include "SpalartAllmaras.H"
makeRASModel(SpalartAllmaras); makeRASModel(SpalartAllmaras);
...@@ -80,12 +87,19 @@ makeRASModel(LRR); ...@@ -80,12 +87,19 @@ makeRASModel(LRR);
makeRASModel(SSG); makeRASModel(SSG);
// -------------------------------------------------------------------------- //
// LES models
// -------------------------------------------------------------------------- //
#include "Smagorinsky.H" #include "Smagorinsky.H"
makeLESModel(Smagorinsky); makeLESModel(Smagorinsky);
#include "WALE.H" #include "WALE.H"
makeLESModel(WALE); makeLESModel(WALE);
#include "dynamicLagrangian.H"
makeLESModel(dynamicLagrangian);
#include "kEqn.H" #include "kEqn.H"
makeLESModel(kEqn); makeLESModel(kEqn);
......
...@@ -61,7 +61,6 @@ SourceFiles ...@@ -61,7 +61,6 @@ SourceFiles
#include "LESeddyViscosity.H" #include "LESeddyViscosity.H"
#include "simpleFilter.H" #include "simpleFilter.H"
#include "LESfilter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-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 "dynamicLagrangian.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace LESModels
{
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class BasicTurbulenceModel>
void dynamicLagrangian<BasicTurbulenceModel>::correctNut
(
const tmp<volTensorField>& gradU
)
{
this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(dev(symm(gradU)));
this->nut_.correctBoundaryConditions();
}
template<class BasicTurbulenceModel>
void dynamicLagrangian<BasicTurbulenceModel>::correctNut()
{
correctNut(fvc::grad(this->U_));
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
dynamicLagrangian<BasicTurbulenceModel>::dynamicLagrangian
(
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
),
flm_
(
IOobject
(
IOobject::groupName("flm", this->U_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
this->mesh_
),
fmm_
(
IOobject
(
IOobject::groupName("fmm", this->U_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
this->mesh_
),
theta_
(
dimensioned<scalar>::lookupOrAddToDict
(
"theta",
this->coeffDict_,
1.5
)
),
simpleFilter_(U.mesh()),
filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
filter_(filterPtr_()),
flm0_("flm0", flm_.dimensions(), 0.0),
fmm0_("fmm0", fmm_.dimensions(), VSMALL)
{
if (type == typeName)
{
correctNut();
this->printCoeffs(type);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
bool dynamicLagrangian<BasicTurbulenceModel>::read()
{
if (LESeddyViscosity<BasicTurbulenceModel>::read())
{
filter_.read(this->coeffDict());
theta_.readIfPresent(this->coeffDict());
return true;
}
else
{
return false;
}
}
template<class BasicTurbulenceModel>
void dynamicLagrangian<BasicTurbulenceModel>::correct()
{
if (!this->turbulence_)
{
return;
}
// Local references
const alphaField& alpha = this->alpha_;
const rhoField& rho = this->rho_;
const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
const volVectorField& U = this->U_;
volScalarField& nut = this->nut_;
LESeddyViscosity<BasicTurbulenceModel>::correct();
tmp<volTensorField> tgradU(fvc::grad(U));
const volTensorField& gradU = tgradU();
volSymmTensorField S(dev(symm(gradU())));
volScalarField magS(mag(S));
volVectorField Uf(filter_(U));
volSymmTensorField Sf(dev(symm(fvc::grad(Uf))));
volScalarField magSf(mag(Sf));
volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
volSymmTensorField M
(
2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
);
volScalarField invT
(
(1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
);
volScalarField LM(L && M);
fvScalarMatrix flmEqn
(
fvm::ddt(flm_)
+ fvm::div(phi(), flm_)
==
invT*LM
- fvm::Sp(invT, flm_)
);
flmEqn.relax();
flmEqn.solve();
bound(flm_, flm0_);
volScalarField MM(M && M);
fvScalarMatrix fmmEqn
(
fvm::ddt(fmm_)
+ fvm::div(phi(), fmm_)
==
invT*MM
- fvm::Sp(invT, fmm_)
);
fmmEqn.relax();
fmmEqn.solve();
bound(fmm_, fmm0_);
correctNut(gradU);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace LESModels
} // End namespace Foam
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-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::dynamicLagrangian
Group
grpLESTurbulence
Description
Dynamic SGS model with Lagrangian averaging
Reference:
\verbatim
Meneveau, C., Lund, T. S., & Cabot, W. H. (1996).
A Lagrangian dynamic subgrid-scale model of turbulence.
Journal of Fluid Mechanics, 319, 353-385.
\endverbatim
SourceFiles
dynamicLagrangian.C
\*---------------------------------------------------------------------------*/
#ifndef dynamicLagrangian_H
#define dynamicLagrangian_H
#include "LESeddyViscosity.H"
#include "simpleFilter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace LESModels
{
/*---------------------------------------------------------------------------*\
Class dynamicLagrangian Declaration
\*---------------------------------------------------------------------------*/
template<class BasicTurbulenceModel>
class dynamicLagrangian
:
public LESeddyViscosity<BasicTurbulenceModel>
{
// Private Member Functions
// Disallow default bitwise copy construct and assignment
dynamicLagrangian(const dynamicLagrangian&);
dynamicLagrangian& operator=(const dynamicLagrangian&);
protected:
// Protected data
volScalarField flm_;
volScalarField fmm_;
dimensionedScalar theta_;
simpleFilter simpleFilter_;
autoPtr<LESfilter> filterPtr_;
LESfilter& filter_;
dimensionedScalar flm0_;
dimensionedScalar fmm0_;
// Protected Member Functions
//- Update sub-grid eddy-viscosity
void correctNut(const tmp<volTensorField>& gradU);
virtual void correctNut();
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("dynamicLagrangian");
// Constructors
//- Construct from components
dynamicLagrangian
(
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 ~dynamicLagrangian()
{}
// Member Functions
//- Read model coefficients if they have changed
virtual bool read();
//- Return SGS kinetic energy
tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
{
return
pow(2.0*flm_/fmm_, 2.0/3.0)
* pow(this->Ce_, -2.0/3.0)
* sqr(this->delta())*magSqr(dev(symm(gradU)));
}
//- Return SGS kinetic energy
virtual tmp<volScalarField> k() const
{
return k(fvc::grad(this->U_));
}
//- Return the effective diffusivity for k
tmp<volScalarField> DkEff() const
{
return tmp<volScalarField>
(
new volScalarField("DkEff", this->nut_ + this->nu())
);
}
//- Correct Eddy-Viscosity and related properties
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace LESModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment