1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 Foam::incompressible::LESModels
28 Namespace for incompressible LES models.
31 Foam::incompressible::LESModel
34 Base class for all incompressible flow LES SGS models.
36 This class defines the basic interface for an incompressible flow SGS
37 model, and encapsulates data of value to all possible models.
38 In particular this includes references to all the dependent fields
39 (U, phi), the physical viscosity nu, and the LESProperties
40 dictionary, which contains the model selection and model coefficients.
45 \*---------------------------------------------------------------------------*/
50 #include "incompressible/turbulenceModel/turbulenceModel.H"
54 #include "fvMatrices.H"
55 #include "incompressible/transportModel/transportModel.H"
56 #include "wallFvPatch.H"
59 #include "runTimeSelectionTables.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 namespace incompressible
68 /*---------------------------------------------------------------------------*\
69 Class LESModel Declaration
70 \*---------------------------------------------------------------------------*/
74 public turbulenceModel,
83 dictionary coeffDict_;
85 dimensionedScalar kMin_;
87 autoPtr<LESdelta> delta_;
90 // Protected Member Functions
92 //- Print model coefficients
93 virtual void printCoeffs();
98 // Private Member Functions
100 //- Disallow default bitwise copy construct
101 LESModel(const LESModel&);
103 //- Disallow default bitwise assignment
104 LESModel& operator=(const LESModel&);
109 //- Runtime type information
110 TypeName("LESModel");
113 // Declare run-time constructor selection table
115 declareRunTimeSelectionTable
121 const volVectorField& U,
122 const surfaceScalarField& phi,
123 transportModel& transport,
124 const word& turbulenceModelName
126 (U, phi, transport, turbulenceModelName)
132 //- Construct from components
136 const volVectorField& U,
137 const surfaceScalarField& phi,
138 transportModel& transport,
139 const word& turbulenceModelName = turbulenceModel::typeName
145 //- Return a reference to the selected LES model
146 static autoPtr<LESModel> New
148 const volVectorField& U,
149 const surfaceScalarField& phi,
150 transportModel& transport,
151 const word& turbulenceModelName = turbulenceModel::typeName
164 //- Const access to the coefficients dictionary,
165 // which provides info. about choice of models,
166 // and all related data (particularly model coefficients).
167 inline const dictionary& coeffDict() const
172 //- Return the lower allowable limit for k (default: SMALL)
173 const dimensionedScalar& kMin() const
178 //- Allow kMin to be changed
179 dimensionedScalar& kMin()
184 //- Access function to filter width
185 virtual const volScalarField& delta() const
191 //- Return the SGS viscosity.
192 virtual tmp<volScalarField> nuSgs() const = 0;
194 //- Return the effective viscosity
195 virtual tmp<volScalarField> nuEff() const
197 return tmp<volScalarField>
199 new volScalarField("nuEff", nuSgs() + nu())
203 //- Return the sub-grid stress tensor.
204 virtual tmp<volSymmTensorField> B() const = 0;
206 //- Return the deviatoric part of the effective sub-grid
207 // turbulence stress tensor including the laminar stress
208 virtual tmp<volSymmTensorField> devBeff() const = 0;
210 //- Returns div(dev(Beff)).
211 // This is the additional term due to the filtering of the NSE.
212 virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const = 0;
215 // RAS compatibility functions for the turbulenceModel base class
217 //- Return the turbulence viscosity
218 virtual tmp<volScalarField> nut() const
223 //- Return the Reynolds stress tensor
224 virtual tmp<volSymmTensorField> R() const
229 //- Return the effective stress tensor including the laminar stress
230 virtual tmp<volSymmTensorField> devReff() const
235 //- Return the source term for the momentum equation
236 virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const
238 return divDevBeff(U);
242 //- Correct Eddy-Viscosity and related properties.
243 // This calls correct(const tmp<volTensorField>& gradU) by supplying
244 // gradU calculated locally.
245 virtual void correct();
247 //- Correct Eddy-Viscosity and related properties
248 virtual void correct(const tmp<volTensorField>& gradU);
250 //- Read LESProperties dictionary
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 } // End namespace incompressible
258 } // End namespace Foam
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 // ************************************************************************* //