1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "dynOneEqEddy.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace compressible
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(dynOneEqEddy, 0);
42 addToRunTimeSelectionTable(LESModel, dynOneEqEddy, dictionary);
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 void dynOneEqEddy::updateSubGridScaleFields(const volSymmTensorField& D)
48 muSgs_ = ck_(D)*rho()*sqrt(k_)*delta();
49 muSgs_.correctBoundaryConditions();
51 alphaSgs_ = muSgs_/Prt_;
52 alphaSgs_.correctBoundaryConditions();
56 dimensionedScalar dynOneEqEddy::ck_(const volSymmTensorField& D) const
58 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
60 volSymmTensorField LL = dev(filter_(sqr(U())) - (sqr(filter_(U()))));
62 volSymmTensorField MM =
63 delta()*(filter_(sqrt(k_)*D) - 2*sqrt(KK + filter_(k_))*filter_(D));
65 return average(LL && MM)/average(magSqr(MM));
69 dimensionedScalar dynOneEqEddy::ce_(const volSymmTensorField& D) const
71 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
74 pow(KK + filter_(k_), 1.5)/(2*delta()) - filter_(pow(k_, 1.5))/delta();
79 filter_(sqrt(k_)*magSqr(D))
80 - 2*sqrt(KK + filter_(k_))*magSqr(filter_(D))
83 return average(ee*mm)/average(mm*mm);
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 dynOneEqEddy::dynOneEqEddy
91 const volScalarField& rho,
92 const volVectorField& U,
93 const surfaceScalarField& phi,
94 const basicThermo& thermoPhysicalModel
97 LESModel(typeName, rho, U, phi, thermoPhysicalModel),
98 GenEddyVisc(rho, U, phi, thermoPhysicalModel),
100 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
101 filter_(filterPtr_())
103 updateSubGridScaleFields(dev(symm(fvc::grad(U))));
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 void dynOneEqEddy::correct(const tmp<volTensorField>& tgradU)
113 const volTensorField& gradU = tgradU();
115 GenEddyVisc::correct(gradU);
117 volSymmTensorField D = dev(symm(gradU));
118 volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
119 volScalarField G = 2*muSgs_*(gradU && D);
124 + fvm::div(phi(), k_)
125 - fvm::laplacian(DkEff(), k_)
128 - fvm::SuSp(2.0/3.0*rho()*divU, k_)
129 - fvm::Sp(ce_(D)*rho()*sqrt(k_)/delta(), k_)
132 bound(k_, dimensionedScalar("0", k_.dimensions(), 1.0e-10));
134 updateSubGridScaleFields(D);
138 bool dynOneEqEddy::read()
140 if (GenEddyVisc::read())
142 filter_.read(coeffDict());
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 } // End namespace LESModels
156 } // End namespace compressible
157 } // End namespace Foam
159 // ************************************************************************* //