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 compressible
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 muSgs_ = ck_(D)*rho()*sqrt(k_)*delta();
48 muSgs_.correctBoundaryConditions();
50 alphaSgs_ = muSgs_/Prt_;
51 alphaSgs_.correctBoundaryConditions();
55 dimensionedScalar homogeneousDynOneEqEddy::ck_(const volSymmTensorField& D) const
57 volScalarField KK(0.5*(filter_(magSqr(U())) - magSqr(filter_(U()))));
59 volSymmTensorField LL(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
63 delta()*(filter_(sqrt(k_)*D) - 2*sqrt(KK + filter_(k_))*filter_(D))
66 return average(LL && MM)/average(magSqr(MM));
70 dimensionedScalar homogeneousDynOneEqEddy::ce_(const volSymmTensorField& D) const
72 volScalarField KK(0.5*(filter_(magSqr(U())) - magSqr(filter_(U()))));
76 pow(KK + filter_(k_), 1.5)/(2*delta()) - filter_(pow(k_, 1.5))/delta()
83 filter_(sqrt(k_)*magSqr(D))
84 - 2*sqrt(KK + filter_(k_))*magSqr(filter_(D))
88 return average(ee*mm)/average(mm*mm);
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 homogeneousDynOneEqEddy::homogeneousDynOneEqEddy
96 const volScalarField& rho,
97 const volVectorField& U,
98 const surfaceScalarField& phi,
99 const basicThermo& thermoPhysicalModel,
100 const word& turbulenceModelName,
101 const word& modelName
104 LESModel(modelName, rho, U, phi, thermoPhysicalModel, turbulenceModelName),
105 GenEddyVisc(rho, U, phi, thermoPhysicalModel),
120 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
121 filter_(filterPtr_())
123 updateSubGridScaleFields(dev(symm(fvc::grad(U))));
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 void homogeneousDynOneEqEddy::correct(const tmp<volTensorField>& tgradU)
133 const volTensorField& gradU = tgradU();
135 GenEddyVisc::correct(gradU);
137 volSymmTensorField D(dev(symm(gradU)));
138 volScalarField divU(fvc::div(phi()/fvc::interpolate(rho())));
139 volScalarField G(2*muSgs_*(gradU && D));
141 tmp<fvScalarMatrix> kEqn
144 + fvm::div(phi(), k_)
145 - fvm::laplacian(DkEff(), k_)
148 - fvm::SuSp(2.0/3.0*rho()*divU, k_)
149 - fvm::Sp(ce_(D)*rho()*sqrt(k_)/delta(), k_)
157 updateSubGridScaleFields(D);
161 bool homogeneousDynOneEqEddy::read()
163 if (GenEddyVisc::read())
165 filter_.read(coeffDict());
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 } // End namespace LESModels
179 } // End namespace compressible
180 } // End namespace Foam
182 // ************************************************************************* //