/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2013-2017 OpenFOAM Foundation
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 .
Namespace
Foam::LESModels
Description
Namespace for LES SGS models.
Class
Foam::LESModel
Description
Templated abstract base class for LES SGS models
SourceFiles
LESModel.C
\*---------------------------------------------------------------------------*/
#ifndef LESModel_H
#define LESModel_H
#include "TurbulenceModel.H"
#include "LESdelta.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class LESModel Declaration
\*---------------------------------------------------------------------------*/
template
class LESModel
:
public BasicTurbulenceModel
{
protected:
// Protected Data
//- LES coefficients dictionary
dictionary LESDict_;
//- Turbulence on/off flag
Switch turbulence_;
//- Flag to print the model coeffs at run-time
Switch printCoeffs_;
//- Model coefficients dictionary
dictionary coeffDict_;
//- Empirical model constant
dimensionedScalar Ce_;
//- Lower limit of k
dimensionedScalar kMin_;
//- Lower limit of epsilon
dimensionedScalar epsilonMin_;
//- Lower limit for omega
dimensionedScalar omegaMin_;
//- Run-time selectable delta model
autoPtr delta_;
// Protected Member Functions
//- Print model coefficients
virtual void printCoeffs(const word& type);
//- No copy construct
LESModel(const LESModel&) = delete;
//- No copy assignment
void operator=(const LESModel&) = delete;
public:
typedef typename BasicTurbulenceModel::alphaField alphaField;
typedef typename BasicTurbulenceModel::rhoField rhoField;
typedef typename BasicTurbulenceModel::transportModel transportModel;
//- Runtime type information
TypeName("LES");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
LESModel,
dictionary,
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
),
(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
);
// Constructors
//- Construct from components
LESModel
(
const word& type,
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
);
// Selectors
//- Return a reference to the selected LES model
static autoPtr New
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName = turbulenceModel::propertiesName
);
//- Destructor
virtual ~LESModel() = default;
// Member Functions
//- Read model coefficients if they have changed
virtual bool read();
// Access
//- Const access to the coefficients dictionary
virtual const dictionary& coeffDict() const
{
return coeffDict_;
}
// Return the empirical model constant
const dimensionedScalar& Ce() const noexcept
{
return Ce_;
}
//- Return the lower allowable limit for k (default: SMALL)
const dimensionedScalar& kMin() const
{
return kMin_;
}
//- Allow kMin to be changed
dimensionedScalar& kMin()
{
return kMin_;
}
//- Access function to filter width
inline const volScalarField& delta() const
{
return *delta_;
}
//- Return the effective viscosity
virtual tmp nuEff() const
{
return tmp
(
new volScalarField
(
IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
this->nut() + this->nu()
)
);
}
//- Return the effective viscosity on patch
virtual tmp nuEff(const label patchi) const
{
return this->nut(patchi) + this->nu(patchi);
}
//- Return the turbulence kinetic energy dissipation rate
virtual tmp epsilon() const;
//- Return the specific dissipation rate
virtual tmp omega() const;
//- Solve the turbulence equations and correct the turbulence viscosity
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "LESModel.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //