1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "dynSmagorinsky.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(dynSmagorinsky, 0);
41 addToRunTimeSelectionTable(LESModel, dynSmagorinsky, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void dynSmagorinsky::updateSubGridScaleFields(const volSymmTensorField& D)
47 nuSgs_ = cD(D)*sqr(delta())*sqrt(magSqr(D));
48 nuSgs_.correctBoundaryConditions();
52 dimensionedScalar dynSmagorinsky::cD(const volSymmTensorField& D) const
54 volSymmTensorField LL = dev(filter_(sqr(U())) - (sqr(filter_(U()))));
56 volSymmTensorField MM =
57 sqr(delta())*(filter_(mag(D)*(D)) - 4*mag(filter_(D))*filter_(D));
59 dimensionedScalar MMMM = average(magSqr(MM));
61 if (MMMM.value() > VSMALL)
63 return average(LL && MM)/MMMM;
72 dimensionedScalar dynSmagorinsky::cI(const volSymmTensorField& D) const
74 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
77 sqr(delta())*(4*sqr(mag(filter_(D))) - filter_(sqr(mag(D))));
79 dimensionedScalar mmmm = average(magSqr(mm));
81 if (mmmm.value() > VSMALL)
83 return average(KK*mm)/mmmm;
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 dynSmagorinsky::dynSmagorinsky
96 const volVectorField& U,
97 const surfaceScalarField& phi,
98 transportModel& transport
101 LESModel(typeName, U, phi, transport),
102 GenEddyVisc(U, phi, transport),
117 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
118 filter_(filterPtr_())
120 updateSubGridScaleFields(dev(symm(fvc::grad(U))));
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 void dynSmagorinsky::correct(const tmp<volTensorField>& gradU)
130 LESModel::correct(gradU);
132 volSymmTensorField D = dev(symm(gradU));
134 k_ = cI(D)*sqr(delta())*magSqr(D);
136 updateSubGridScaleFields(D);
140 bool dynSmagorinsky::read()
142 if (GenEddyVisc::read())
144 filter_.read(coeffDict());
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 } // End namespace LESModels
158 } // End namespace incompressible
159 } // End namespace Foam
161 // ************************************************************************* //