initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / incompressible / locDynOneEqEddy / locDynOneEqEddy.C
blob35f99986add4090ce9b6cca11716ac528bca8a6e
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 "locDynOneEqEddy.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(locDynOneEqEddy, 0);
42 addToRunTimeSelectionTable(LESModel, locDynOneEqEddy, dictionary);
44 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
46 volScalarField locDynOneEqEddy::ck
48     const volSymmTensorField& D,
49     const volScalarField& KK
50 ) const
52     volSymmTensorField LL =
53         simpleFilter_(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
55     volSymmTensorField MM = simpleFilter_(-2.0*delta()*pow(KK, 0.5)*filter_(D));
57     volScalarField ck =
58         simpleFilter_(0.5*(LL && MM))
59        /(
60             simpleFilter_(magSqr(MM))
61           + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
62         );
64     return 0.5*(mag(ck) + ck);
68 volScalarField locDynOneEqEddy::ce
70     const volSymmTensorField& D,
71     const volScalarField& KK
72 ) const
74     volScalarField ce =
75         simpleFilter_(nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
76        /simpleFilter_(pow(KK, 1.5)/(2.0*delta()));
78     return 0.5*(mag(ce) + ce);
82 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
84 locDynOneEqEddy::locDynOneEqEddy
86     const volVectorField& U,
87     const surfaceScalarField& phi,
88     transportModel& transport
91     LESModel(typeName, U, phi, transport),
92     GenEddyVisc(U, phi, transport),
94     k_
95     (
96         IOobject
97         (
98             "k",
99             runTime_.timeName(),
100             mesh_,
101             IOobject::MUST_READ,
102             IOobject::AUTO_WRITE
103         ),
104         mesh_
105     ),
107     simpleFilter_(U.mesh()),
108     filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
109     filter_(filterPtr_())
111     printCoeffs();
115 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
117 locDynOneEqEddy::~locDynOneEqEddy()
121 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
123 void locDynOneEqEddy::correct(const tmp<volTensorField>& gradU)
125     LESModel::correct(gradU);
127     volSymmTensorField D = symm(gradU);
129     volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
130     KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
132     volScalarField P = 2.0*nuSgs_*magSqr(D);
134     solve
135     (
136        fvm::ddt(k_)
137      + fvm::div(phi(), k_)
138      - fvm::laplacian(DkEff(), k_)
139     ==
140        P
141      - fvm::Sp(ce(D, KK)*sqrt(k_)/delta(), k_)
142     );
144     bound(k_, k0());
146     nuSgs_ = ck(D, KK)*sqrt(k_)*delta();
147     nuSgs_.correctBoundaryConditions();
151 bool locDynOneEqEddy::read()
153     if (GenEddyVisc::read())
154     {
155         filter_.read(coeffDict());
157         return true;
158     }
159     else
160     {
161         return false;
162     }
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 } // End namespace LESModels
169 } // End namespace incompressible
170 } // End namespace Foam
172 // ************************************************************************* //