correction to d4edb38234db8268907f04836d49bb93461b8a88
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / incompressible / Smagorinsky2 / Smagorinsky2.C
blob351db3526e7d869f6df6dfd54d9e712b7021c3c5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "Smagorinsky2.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(Smagorinsky2, 0);
42 addToRunTimeSelectionTable(LESModel, Smagorinsky2, dictionary);
44 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
46 Smagorinsky2::Smagorinsky2
48     const volVectorField& U,
49     const surfaceScalarField& phi,
50     transportModel& transport
53     LESModel(typeName, U, phi, transport),
54     Smagorinsky(U, phi, transport),
56     cD2_
57     (
58         dimensioned<scalar>::lookupOrAddToDict
59         (
60             "cD2",
61             coeffDict(),
62             0.02
63         )
64     )
66     printCoeffs();
70 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
72 // Evaluate B (from the definition of an eddy viscosity model) and
73 // return it.
75 tmp<volSymmTensorField> Smagorinsky2::B() const
77     volSymmTensorField D = dev(symm(fvc::grad(U())));
79     return (((2.0/3.0)*I)*k() - 2.0*nuSgs_*D - (2.0*cD2_)*delta()*(D&D));
83 tmp<fvVectorMatrix> Smagorinsky2::divDevBeff
85     volVectorField& U
86 ) const
88     volTensorField gradU = fvc::grad(U);
90     volSymmTensorField aniNuEff
91     (
92         "aniNuEff",
93         I*nuEff() + cD2_*delta()*symm(gradU)
94     );
96     return
97     (
98       - fvm::laplacian(aniNuEff, U) - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
99     );
103 bool Smagorinsky2::read()
105     if (Smagorinsky::read())
106     {
107         cD2.readIfPresent(coeffDict());
109         return true;
110     }
111     else
112     {
113         return false;
114     }
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 } // End namespace LESModels
121 } // End namespace incompressible
122 } // End namespace Foam
124 // ************************************************************************* //