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 "homogeneousDynOneEqEddy.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(homogeneousDynOneEqEddy, 0);
41 addToRunTimeSelectionTable(LESModel, homogeneousDynOneEqEddy, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void homogeneousDynOneEqEddy::updateSubGridScaleFields(const volSymmTensorField& D)
47 nuSgs_ = ck(D)*sqrt(k_)*delta();
48 nuSgs_.correctBoundaryConditions();
52 dimensionedScalar homogeneousDynOneEqEddy::ck(const volSymmTensorField& D) const
54 tmp<volScalarField> KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
56 const volSymmTensorField MM
58 delta()*(filter_(sqrt(k_)*D) - 2*sqrt(KK + filter_(k_))*filter_(D))
61 dimensionedScalar MMMM = average(magSqr(MM));
63 if (MMMM.value() > VSMALL)
65 tmp<volSymmTensorField> LL = dev(filter_(sqr(U())) - sqr(filter_(U())));
67 return average(LL && MM)/MMMM;
76 dimensionedScalar homogeneousDynOneEqEddy::ce(const volSymmTensorField& D) const
78 const volScalarField KK
80 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())))
83 const volScalarField mm
85 pow(KK + filter_(k_), 1.5)/(2*delta()) - filter_(pow(k_, 1.5))/delta()
88 dimensionedScalar mmmm = average(magSqr(mm));
90 if (mmmm.value() > VSMALL)
92 tmp<volScalarField> ee =
96 filter_(sqrt(k_)*magSqr(D))
97 - 2*sqrt(KK + filter_(k_))*magSqr(filter_(D))
101 return average(ee*mm)/mmmm;
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 homogeneousDynOneEqEddy::homogeneousDynOneEqEddy
114 const volVectorField& U,
115 const surfaceScalarField& phi,
116 transportModel& transport,
117 const word& turbulenceModelName,
118 const word& modelName
121 LESModel(modelName, U, phi, transport, turbulenceModelName),
122 GenEddyVisc(U, phi, transport),
137 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
138 filter_(filterPtr_())
142 updateSubGridScaleFields(symm(fvc::grad(U)));
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 void homogeneousDynOneEqEddy::correct(const tmp<volTensorField>& tgradU)
152 const volTensorField& gradU = tgradU();
154 GenEddyVisc::correct(gradU);
156 const volSymmTensorField D(symm(gradU));
158 const volScalarField P(2.0*nuSgs_*magSqr(D));
160 tmp<fvScalarMatrix> kEqn
163 + fvm::div(phi(), k_)
164 - fvm::laplacian(DkEff(), k_)
167 - fvm::Sp(ce(D)*sqrt(k_)/delta(), k_)
175 updateSubGridScaleFields(D);
179 bool homogeneousDynOneEqEddy::read()
181 if (GenEddyVisc::read())
183 filter_.read(coeffDict());
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 } // End namespace LESModels
197 } // End namespace incompressible
198 } // End namespace Foam
200 // ************************************************************************* //