Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / Smagorinsky2 / Smagorinsky2.C
blob97d32bce9dd2ac39bc2fb3184487085be42132e9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "Smagorinsky2.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace incompressible
35 namespace LESModels
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(Smagorinsky2, 0);
41 addToRunTimeSelectionTable(LESModel, Smagorinsky2, dictionary);
43 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
45 Smagorinsky2::Smagorinsky2
47     const volVectorField& U,
48     const surfaceScalarField& phi,
49     transportModel& transport,
50     const word& turbulenceModelName,
51     const word& modelName
54     LESModel(modelName, U, phi, transport, turbulenceModelName),
55     Smagorinsky(U, phi, transport),
57     cD2_
58     (
59         dimensioned<scalar>::lookupOrAddToDict
60         (
61             "cD2",
62             coeffDict_,
63             0.02
64         )
65     )
67     printCoeffs();
71 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
73 // Evaluate B (from the definition of an eddy viscosity model) and
74 // return it.
76 tmp<volSymmTensorField> Smagorinsky2::B() const
78     volSymmTensorField D(dev(symm(fvc::grad(U()))));
80     return (((2.0/3.0)*I)*k() - 2.0*nuSgs_*D - (2.0*cD2_)*delta()*(D&D));
84 tmp<fvVectorMatrix> Smagorinsky2::divDevBeff
86     volVectorField& U
87 ) const
89     volTensorField gradU(fvc::grad(U));
91     volSymmTensorField aniNuEff
92     (
93         "aniNuEff",
94         I*nuEff() + cD2_*delta()*symm(gradU)
95     );
97     return
98     (
99       - fvm::laplacian(aniNuEff, U) - fvc::div(nuEff()*dev(T(fvc::grad(U))))
100     );
104 bool Smagorinsky2::read()
106     if (Smagorinsky::read())
107     {
108         cD2.readIfPresent(coeffDict());
110         return true;
111     }
112     else
113     {
114         return false;
115     }
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 } // End namespace LESModels
122 } // End namespace incompressible
123 } // End namespace Foam
125 // ************************************************************************* //