initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / dynSmagorinsky / dynSmagorinsky.C
blobfedd18ee7d169d5fcb7f2741dc1e05f82f78db2a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 "dynSmagorinsky.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(dynSmagorinsky, 0);
42 addToRunTimeSelectionTable(LESModel, dynSmagorinsky, dictionary);
44 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
46 void dynSmagorinsky::updateSubGridScaleFields(const volSymmTensorField& D)
48     nuSgs_ = cD(D)*sqr(delta())*sqrt(magSqr(D));
49     nuSgs_.correctBoundaryConditions();
53 dimensionedScalar dynSmagorinsky::cD(const volSymmTensorField& D) const
55     volSymmTensorField LL = dev(filter_(sqr(U())) - (sqr(filter_(U()))));
57     volSymmTensorField MM =
58         sqr(delta())*(filter_(mag(D)*(D)) - 4*mag(filter_(D))*filter_(D));
60     dimensionedScalar MMMM = average(magSqr(MM));
62     if (MMMM.value() > VSMALL)
63     {
64         return average(LL && MM)/MMMM;
65     }
66     else
67     {
68         return 0.0;
69     }
73 dimensionedScalar dynSmagorinsky::cI(const volSymmTensorField& D) const
75     volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
77     volScalarField mm =
78         sqr(delta())*(4*sqr(mag(filter_(D))) - filter_(sqr(mag(D))));
80     dimensionedScalar mmmm = average(magSqr(mm));
82     if (mmmm.value() > VSMALL)
83     {
84         return average(KK*mm)/mmmm;
85     }
86     else
87     {
88         return 0.0;
89     }
93 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
95 dynSmagorinsky::dynSmagorinsky
97     const volVectorField& U,
98     const surfaceScalarField& phi,
99     transportModel& transport
102     LESModel(typeName, U, phi, transport),
103     GenEddyVisc(U, phi, transport),
105     k_
106     (
107         IOobject
108         (
109             "k",
110             runTime_.timeName(),
111             mesh_,
112             IOobject::MUST_READ,
113             IOobject::AUTO_WRITE
114         ),
115         mesh_
116     ),
118     filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
119     filter_(filterPtr_())
121     updateSubGridScaleFields(dev(symm(fvc::grad(U))));
123     printCoeffs();
127 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
129 void dynSmagorinsky::correct(const tmp<volTensorField>& gradU)
131     LESModel::correct(gradU);
133     volSymmTensorField D = dev(symm(gradU));
135     k_ = cI(D)*sqr(delta())*magSqr(D);
137     updateSubGridScaleFields(D);
141 bool dynSmagorinsky::read()
143     if (GenEddyVisc::read())
144     {
145         filter_.read(coeffDict());
147         return true;
148     }
149     else
150     {
151         return false;
152     }
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 } // End namespace LESModels
159 } // End namespace incompressible
160 } // End namespace Foam
162 // ************************************************************************* //