Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / GenEddyVisc / GenEddyVisc.C
blobbe9eb9b228bc2f651febb00d29cb497a8ec3a0c2
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 "GenEddyVisc.H"
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 namespace Foam
32 namespace incompressible
34 namespace LESModels
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 GenEddyVisc::GenEddyVisc
41     const volVectorField& U,
42     const surfaceScalarField& phi,
43     transportModel& transport
46     LESModel(word("GenEddyVisc"), U, phi, transport),
48     ce_
49     (
50         dimensionedScalar::lookupOrAddToDict
51         (
52             "ce",
53             coeffDict_,
54             1.048
55         )
56     ),
58     nuSgs_
59     (
60         IOobject
61         (
62             "nuSgs",
63             runTime_.timeName(),
64             U_.db(),
65             IOobject::MUST_READ,
66             IOobject::AUTO_WRITE
67         ),
68         mesh_
69     )
71 //    printCoeffs();
75 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
77 tmp<volSymmTensorField> GenEddyVisc::B() const
79     return ((2.0/3.0)*I)*k() - nuSgs_*twoSymm(fvc::grad(U()));
83 tmp<volSymmTensorField> GenEddyVisc::devBeff() const
85     return -nuEff()*dev(twoSymm(fvc::grad(U())));
89 tmp<fvVectorMatrix> GenEddyVisc::divDevBeff(volVectorField& U) const
91     return
92     (
93       - fvm::laplacian(nuEff(), U) - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
94     );
98 void GenEddyVisc::correct(const tmp<volTensorField>& gradU)
100     LESModel::correct(gradU);
104 bool GenEddyVisc::read()
106     if (LESModel::read())
107     {
108         ce_.readIfPresent(coeffDict());
110         return true;
111     }
112     else
113     {
114         return false;
115     }
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 } // End namespace LESModels
122 } // End namespace incompressible
123 } // End namespace Foam
125 // ************************************************************************* //