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::oneEqEddy
28 One Equation Eddy Viscosity Model for incompressible flows
30 Eddy viscosity SGS model using a modeled balance equation to simulate the
31 behaviour of k, hence,
33 d/dt(k) + div(U*k) - div(nuEff*grad(k))
39 B = 2/3*k*I - 2*nuSgs*dev(D)
40 Beff = 2/3*k*I - 2*nuEff*dev(D)
45 nuSgs = ck*sqrt(k)*delta
52 \*---------------------------------------------------------------------------*/
57 #include "GenEddyVisc.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 namespace incompressible
68 /*---------------------------------------------------------------------------*\
69 Class oneEqEddy Declaration
70 \*---------------------------------------------------------------------------*/
80 dimensionedScalar ck_;
83 // Private Member Functions
85 //- Update sub-grid scale fields
86 void updateSubGridScaleFields();
88 // Disallow default bitwise copy construct and assignment
89 oneEqEddy(const oneEqEddy&);
90 oneEqEddy& operator=(const oneEqEddy&);
95 //- Runtime type information
96 TypeName("oneEqEddy");
100 //- Construct from components
103 const volVectorField& U,
104 const surfaceScalarField& phi,
105 transportModel& transport,
106 const word& turbulenceModelName = turbulenceModel::typeName,
107 const word& modelName = typeName
118 //- Return SGS kinetic energy
119 virtual tmp<volScalarField> k() const
124 //- Return the effective diffusivity for k
125 tmp<volScalarField> DkEff() const
127 return tmp<volScalarField>
129 new volScalarField("DkEff", nuSgs_ + nu())
133 //- Correct Eddy-Viscosity and related properties
134 virtual void correct(const tmp<volTensorField>& gradU);
136 //- Read LESProperties dictionary
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 } // End namespace LESModels
144 } // End namespace incompressible
145 } // End namespace Foam
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 // ************************************************************************* //