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 "locDynOneEqEddy.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(locDynOneEqEddy, 0);
41 addToRunTimeSelectionTable(LESModel, locDynOneEqEddy, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void locDynOneEqEddy::updateSubGridScaleFields
47 const volSymmTensorField& D,
48 const volScalarField& KK
51 nuSgs_ = ck(D, KK)*sqrt(k_)*delta();
52 nuSgs_.correctBoundaryConditions();
56 volScalarField locDynOneEqEddy::ck
58 const volSymmTensorField& D,
59 const volScalarField& KK
62 volSymmTensorField LL =
63 simpleFilter_(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
65 volSymmTensorField MM = simpleFilter_(-2.0*delta()*pow(KK, 0.5)*filter_(D));
68 simpleFilter_(0.5*(LL && MM))
70 simpleFilter_(magSqr(MM))
71 + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
74 return 0.5*(mag(ck) + ck);
78 volScalarField locDynOneEqEddy::ce
80 const volSymmTensorField& D,
81 const volScalarField& KK
85 simpleFilter_(nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
86 /simpleFilter_(pow(KK, 1.5)/(2.0*delta()));
88 return 0.5*(mag(ce) + ce);
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 locDynOneEqEddy::locDynOneEqEddy
96 const volVectorField& U,
97 const surfaceScalarField& phi,
98 transportModel& transport
101 LESModel(typeName, U, phi, transport),
102 GenEddyVisc(U, phi, transport),
117 simpleFilter_(U.mesh()),
118 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
119 filter_(filterPtr_())
121 volScalarField KK = 0.5*(filter_(magSqr(U)) - magSqr(filter_(U)));
122 updateSubGridScaleFields(symm(fvc::grad(U)), KK);
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 void locDynOneEqEddy::correct(const tmp<volTensorField>& gradU)
132 LESModel::correct(gradU);
134 volSymmTensorField D = symm(gradU);
136 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
137 KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
139 volScalarField P = 2.0*nuSgs_*magSqr(D);
144 + fvm::div(phi(), k_)
145 - fvm::laplacian(DkEff(), k_)
148 - fvm::Sp(ce(D, KK)*sqrt(k_)/delta(), k_)
156 updateSubGridScaleFields(D, KK);
160 bool locDynOneEqEddy::read()
162 if (GenEddyVisc::read())
164 filter_.read(coeffDict());
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 } // End namespace LESModels
178 } // End namespace incompressible
179 } // End namespace Foam
181 // ************************************************************************* //