initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / incompressible / dynMixedSmagorinsky / dynMixedSmagorinsky.C
blobaabecafd199797542a6e2ccfa1d7449b314ef0b9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 "dynMixedSmagorinsky.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(dynMixedSmagorinsky, 0);
42 addToRunTimeSelectionTable(LESModel, dynMixedSmagorinsky, dictionary);
44 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
46 dynMixedSmagorinsky::dynMixedSmagorinsky
48     const volVectorField& U,
49     const surfaceScalarField& phi,
50     transportModel& transport
53     LESModel(typeName, U, phi, transport),
54     scaleSimilarity(U, phi, transport),
55     dynSmagorinsky(U, phi, transport)
57     printCoeffs();
61 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
63 tmp<volScalarField> dynMixedSmagorinsky::k() const
65     return
66     (
67         scaleSimilarity::k()
68       + dynSmagorinsky::k()
69     );
73 tmp<volScalarField> dynMixedSmagorinsky::epsilon() const
75     return
76     (
77         scaleSimilarity::epsilon()
78       + dynSmagorinsky::epsilon()
79     );
83 tmp<volSymmTensorField> dynMixedSmagorinsky::B() const
85     return
86     (
87         scaleSimilarity::B()
88       + dynSmagorinsky::B()
89     );
93 tmp<volSymmTensorField> dynMixedSmagorinsky::devBeff() const
95     return
96     (
97         scaleSimilarity::devBeff()
98       + dynSmagorinsky::devBeff()
99     );
103 tmp<fvVectorMatrix> dynMixedSmagorinsky::divDevBeff(volVectorField& U) const
105     return
106     (
107         scaleSimilarity::divDevBeff(U)
108       + dynSmagorinsky::divDevBeff(U)
109     );
113 void dynMixedSmagorinsky::correct(const tmp<volTensorField>& gradU)
115     scaleSimilarity::correct(gradU);
116     dynSmagorinsky::correct(gradU);
120 bool dynMixedSmagorinsky::read()
122     if (LESModel::read())
123     {
124         scaleSimilarity::read();
125         dynSmagorinsky::read();
127         return true;
128     }
129     else
130     {
131         return false;
132     }
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 } // End namespace LESModels
139 } // namespace incompressible
140 } // End namespace Foam
142 // ************************************************************************* //