Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / spectEddyVisc / spectEddyVisc.C
blobce562c40767f3fbe7f447e2e97be3f3c1ebe5e49
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace incompressible
35 namespace LESModels
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++)
49     {
50         nuSgs_ =
51             nu()
52            /(
53                  scalar(1)
54                - exp(-cB_*pow(nu()/(nuSgs_ + nu()), 1.0/3.0)*pow(Re, -2.0/3.0))
55             );
56     }
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),
75     cB_
76     (
77         dimensionedScalar::lookupOrAddToDict
78         (
79             "cB",
80             coeffDict_,
81             8.22
82         )
83     ),
84     cK1_
85     (
86         dimensionedScalar::lookupOrAddToDict
87         (
88             "cK1",
89             coeffDict_,
90             0.83
91         )
92     ),
93     cK2_
94     (
95         dimensionedScalar::lookupOrAddToDict
96         (
97             "cK2",
98             coeffDict_,
99             1.03
100         )
101     ),
102     cK3_
103     (
104         dimensionedScalar::lookupOrAddToDict
105         (
106             "cK3",
107             coeffDict_,
108             4.75
109         )
110     ),
111     cK4_
112     (
113         dimensionedScalar::lookupOrAddToDict
114         (
115             "cK4",
116             coeffDict_,
117             2.55
118         )
119     )
121     printCoeffs();
123     updateSubGridScaleFields(fvc::grad(U));
127 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
129 tmp<volScalarField> spectEddyVisc::k() const
131     volScalarField eps = 2*nuEff()*magSqr(symm(fvc::grad(U())));
133     return
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())
151     {
152         cB_.readIfPresent(coeffDict());
153         cK1_.readIfPresent(coeffDict());
154         cK2_.readIfPresent(coeffDict());
155         cK3_.readIfPresent(coeffDict());
156         cK4_.readIfPresent(coeffDict());
158         return true;
159     }
160     else
161     {
162         return false;
163     }
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 } // End namespace LESModels
170 } // End namespace incompressible
171 } // End namespace Foam
173 // ************************************************************************* //