initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / locDynOneEqEddy / locDynOneEqEddy.C
blobf463e3f2082aad874af2ce36dfd8def76e93ed75
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 void locDynOneEqEddy::updateSubGridScaleFields
48     const volSymmTensorField& D,
49     const volScalarField& KK
52     nuSgs_ = ck(D, KK)*sqrt(k_)*delta();
53     nuSgs_.correctBoundaryConditions();
57 volScalarField locDynOneEqEddy::ck
59     const volSymmTensorField& D,
60     const volScalarField& KK
61 ) const
63     volSymmTensorField LL =
64         simpleFilter_(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
66     volSymmTensorField MM = simpleFilter_(-2.0*delta()*pow(KK, 0.5)*filter_(D));
68     volScalarField ck =
69         simpleFilter_(0.5*(LL && MM))
70        /(
71             simpleFilter_(magSqr(MM))
72           + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
73         );
75     return 0.5*(mag(ck) + ck);
79 volScalarField locDynOneEqEddy::ce
81     const volSymmTensorField& D,
82     const volScalarField& KK
83 ) const
85     volScalarField ce =
86         simpleFilter_(nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
87        /simpleFilter_(pow(KK, 1.5)/(2.0*delta()));
89     return 0.5*(mag(ce) + ce);
93 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
95 locDynOneEqEddy::locDynOneEqEddy
97     const volVectorField& U,
98     const surfaceScalarField& phi,
99     transportModel& transport
102     LESModel(typeName, U, phi, transport),
103     GenEddyVisc(U, phi, transport),
105     k_
106     (
107         IOobject
108         (
109             "k",
110             runTime_.timeName(),
111             mesh_,
112             IOobject::MUST_READ,
113             IOobject::AUTO_WRITE
114         ),
115         mesh_
116     ),
118     simpleFilter_(U.mesh()),
119     filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
120     filter_(filterPtr_())
122     volScalarField KK = 0.5*(filter_(magSqr(U)) - magSqr(filter_(U)));
123     updateSubGridScaleFields(symm(fvc::grad(U)), KK);
125     printCoeffs();
129 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
131 void locDynOneEqEddy::correct(const tmp<volTensorField>& gradU)
133     LESModel::correct(gradU);
135     volSymmTensorField D = symm(gradU);
137     volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
138     KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
140     volScalarField P = 2.0*nuSgs_*magSqr(D);
142     fvScalarMatrix kEqn
143     (
144        fvm::ddt(k_)
145      + fvm::div(phi(), k_)
146      - fvm::laplacian(DkEff(), k_)
147     ==
148        P
149      - fvm::Sp(ce(D, KK)*sqrt(k_)/delta(), k_)
150     );
152     kEqn.relax();
153     kEqn.solve();
155     bound(k_, k0());
157     updateSubGridScaleFields(D, KK);
161 bool locDynOneEqEddy::read()
163     if (GenEddyVisc::read())
164     {
165         filter_.read(coeffDict());
167         return true;
168     }
169     else
170     {
171         return false;
172     }
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 } // End namespace LESModels
179 } // End namespace incompressible
180 } // End namespace Foam
182 // ************************************************************************* //