Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / dynMixedSmagorinsky / dynMixedSmagorinsky.C
blobbdc1e7e0eef5d56728fb20b0275186149a043635
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 "dynMixedSmagorinsky.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace incompressible
35 namespace LESModels
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(dynMixedSmagorinsky, 0);
41 addToRunTimeSelectionTable(LESModel, dynMixedSmagorinsky, dictionary);
43 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
45 dynMixedSmagorinsky::dynMixedSmagorinsky
47     const volVectorField& U,
48     const surfaceScalarField& phi,
49     transportModel& transport
52     LESModel(typeName, U, phi, transport),
53     scaleSimilarity(U, phi, transport),
54     dynSmagorinsky(U, phi, transport)
56     printCoeffs();
60 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
62 tmp<volScalarField> dynMixedSmagorinsky::k() const
64     return
65     (
66         scaleSimilarity::k()
67       + dynSmagorinsky::k()
68     );
72 tmp<volScalarField> dynMixedSmagorinsky::epsilon() const
74     return
75     (
76         scaleSimilarity::epsilon()
77       + dynSmagorinsky::epsilon()
78     );
82 tmp<volSymmTensorField> dynMixedSmagorinsky::B() const
84     return
85     (
86         scaleSimilarity::B()
87       + dynSmagorinsky::B()
88     );
92 tmp<volSymmTensorField> dynMixedSmagorinsky::devBeff() const
94     return
95     (
96         scaleSimilarity::devBeff()
97       + dynSmagorinsky::devBeff()
98     );
102 tmp<fvVectorMatrix> dynMixedSmagorinsky::divDevBeff(volVectorField& U) const
104     return
105     (
106         scaleSimilarity::divDevBeff(U)
107       + dynSmagorinsky::divDevBeff(U)
108     );
112 void dynMixedSmagorinsky::correct(const tmp<volTensorField>& gradU)
114     scaleSimilarity::correct(gradU);
115     dynSmagorinsky::correct(gradU());
119 bool dynMixedSmagorinsky::read()
121     if (LESModel::read())
122     {
123         scaleSimilarity::read();
124         dynSmagorinsky::read();
126         return true;
127     }
128     else
129     {
130         return false;
131     }
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 } // End namespace LESModels
138 } // namespace incompressible
139 } // End namespace Foam
141 // ************************************************************************* //