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 "spectEddyVisc.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(spectEddyVisc, 0);
41 addToRunTimeSelectionTable(LESModel, spectEddyVisc, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void spectEddyVisc::updateSubGridScaleFields(const volTensorField& gradU)
47 volScalarField Re = sqr(delta())*mag(symm(gradU))/nu();
48 for (label i=0; i<5; i++)
54 - exp(-cB_*pow(nu()/(nuSgs_ + nu()), 1.0/3.0)*pow(Re, -2.0/3.0))
58 nuSgs_.correctBoundaryConditions();
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 spectEddyVisc::spectEddyVisc
66 const volVectorField& U,
67 const surfaceScalarField& phi,
68 transportModel& transport
71 LESModel(typeName, U, phi, transport),
72 GenEddyVisc(U, phi, transport),
77 dimensionedScalar::lookupOrAddToDict
86 dimensionedScalar::lookupOrAddToDict
95 dimensionedScalar::lookupOrAddToDict
104 dimensionedScalar::lookupOrAddToDict
113 dimensionedScalar::lookupOrAddToDict
123 updateSubGridScaleFields(fvc::grad(U));
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 tmp<volScalarField> spectEddyVisc::k() const
131 volScalarField eps = 2*nuEff()*magSqr(symm(fvc::grad(U())));
134 cK1_*pow(delta()*eps, 2.0/3.0)
135 *exp(-cK2_*pow(delta(), -4.0/3.0)*nu()/pow(eps, 1.0/3.0))
136 - cK3_*sqrt(eps*nu())
137 *erfc(cK4_*pow(delta(), -2.0/3.0)*sqrt(nu())*pow(eps, -1.0/6.0));
141 void spectEddyVisc::correct(const tmp<volTensorField>& gradU)
143 GenEddyVisc::correct(gradU);
144 updateSubGridScaleFields(gradU());
148 bool spectEddyVisc::read()
150 if (GenEddyVisc::read())
152 cB_.readIfPresent(coeffDict());
153 cK1_.readIfPresent(coeffDict());
154 cK2_.readIfPresent(coeffDict());
155 cK3_.readIfPresent(coeffDict());
156 cK4_.readIfPresent(coeffDict());
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 } // End namespace LESModels
170 } // End namespace incompressible
171 } // End namespace Foam
173 // ************************************************************************* //