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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "homogeneousDynSmagorinsky.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(homogeneousDynSmagorinsky, 0);
41 addToRunTimeSelectionTable(LESModel, homogeneousDynSmagorinsky, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void homogeneousDynSmagorinsky::updateSubGridScaleFields
47 const volSymmTensorField& D
50 nuSgs_ = cD(D)*sqr(delta())*sqrt(magSqr(D));
51 nuSgs_.correctBoundaryConditions();
55 dimensionedScalar homogeneousDynSmagorinsky::cD
57 const volSymmTensorField& D
60 const volSymmTensorField MM
62 sqr(delta())*(filter_(mag(D)*(D)) - 4*mag(filter_(D))*filter_(D))
65 dimensionedScalar MMMM = average(magSqr(MM));
67 if (MMMM.value() > VSMALL)
69 tmp<volSymmTensorField> LL =
70 dev(filter_(sqr(U())) - (sqr(filter_(U()))));
72 return 0.5*average(LL && MM)/MMMM;
81 dimensionedScalar homogeneousDynSmagorinsky::cI
83 const volSymmTensorField& D
86 const volScalarField mm
88 sqr(delta())*(4*sqr(mag(filter_(D))) - filter_(sqr(mag(D))))
91 dimensionedScalar mmmm = average(magSqr(mm));
93 if (mmmm.value() > VSMALL)
95 tmp<volScalarField> KK =
96 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
98 return average(KK*mm)/mmmm;
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 homogeneousDynSmagorinsky::homogeneousDynSmagorinsky
111 const volVectorField& U,
112 const surfaceScalarField& phi,
113 transportModel& transport,
114 const word& turbulenceModelName,
115 const word& modelName
118 LESModel(modelName, U, phi, transport, turbulenceModelName),
119 GenEddyVisc(U, phi, transport),
134 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
135 filter_(filterPtr_())
139 updateSubGridScaleFields(dev(symm(fvc::grad(U))));
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 void homogeneousDynSmagorinsky::correct(const tmp<volTensorField>& gradU)
149 LESModel::correct(gradU);
151 const volSymmTensorField D(dev(symm(gradU)));
153 k_ = cI(D)*sqr(delta())*magSqr(D);
156 updateSubGridScaleFields(D);
160 bool homogeneousDynSmagorinsky::read()
162 if (GenEddyVisc::read())
164 filter_.read(coeffDict());
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 } // End namespace LESModels
178 } // End namespace incompressible
179 } // End namespace Foam
181 // ************************************************************************* //