FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

dynOneEqEddy.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "dynOneEqEddy.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace incompressible
00034 {
00035 namespace LESModels
00036 {
00037 
00038 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(dynOneEqEddy, 0);
00041 addToRunTimeSelectionTable(LESModel, dynOneEqEddy, dictionary);
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void dynOneEqEddy::updateSubGridScaleFields(const volSymmTensorField& D)
00046 {
00047     nuSgs_ = ck(D)*sqrt(k_)*delta();
00048     nuSgs_.correctBoundaryConditions();
00049 }
00050 
00051 
00052 dimensionedScalar dynOneEqEddy::ck(const volSymmTensorField& D) const
00053 {
00054     volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
00055 
00056     volSymmTensorField LL = dev(filter_(sqr(U())) - sqr(filter_(U())));
00057 
00058     volSymmTensorField MM =
00059         delta()*(filter_(sqrt(k_)*D) - 2*sqrt(KK + filter_(k_))*filter_(D));
00060 
00061     dimensionedScalar MMMM = average(magSqr(MM));
00062 
00063     if (MMMM.value() > VSMALL)
00064     {
00065         return average(LL && MM)/MMMM;
00066     }
00067     else
00068     {
00069         return 0.0;
00070     }
00071 }
00072 
00073 
00074 dimensionedScalar dynOneEqEddy::ce(const volSymmTensorField& D) const
00075 {
00076     volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
00077 
00078     volScalarField mm =
00079         pow(KK + filter_(k_), 1.5)/(2*delta()) - filter_(pow(k_, 1.5))/delta();
00080 
00081     volScalarField ee =
00082         2*delta()*ck(D)
00083        *(
00084            filter_(sqrt(k_)*magSqr(D))
00085          - 2*sqrt(KK + filter_(k_))*magSqr(filter_(D))
00086         );
00087 
00088     dimensionedScalar mmmm = average(magSqr(mm));
00089 
00090     if (mmmm.value() > VSMALL)
00091     {
00092         return average(ee*mm)/mmmm;
00093     }
00094     else
00095     {
00096         return 0.0;
00097     }
00098 }
00099 
00100 
00101 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00102 
00103 dynOneEqEddy::dynOneEqEddy
00104 (
00105     const volVectorField& U,
00106     const surfaceScalarField& phi,
00107     transportModel& transport
00108 )
00109 :
00110     LESModel(typeName, U, phi, transport),
00111     GenEddyVisc(U, phi, transport),
00112 
00113     k_
00114     (
00115         IOobject
00116         (
00117             "k",
00118             runTime_.timeName(),
00119             mesh_,
00120             IOobject::MUST_READ,
00121             IOobject::AUTO_WRITE
00122         ),
00123         mesh_
00124     ),
00125 
00126     filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
00127     filter_(filterPtr_())
00128 {
00129     updateSubGridScaleFields(symm(fvc::grad(U)));
00130 
00131     printCoeffs();
00132 }
00133 
00134 
00135 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00136 
00137 void dynOneEqEddy::correct(const tmp<volTensorField>& gradU)
00138 {
00139     GenEddyVisc::correct(gradU);
00140 
00141     volSymmTensorField D = symm(gradU);
00142 
00143     volScalarField P = 2.0*nuSgs_*magSqr(D);
00144 
00145     fvScalarMatrix kEqn
00146     (
00147        fvm::ddt(k_)
00148      + fvm::div(phi(), k_)
00149      - fvm::laplacian(DkEff(), k_)
00150     ==
00151        P
00152      - fvm::Sp(ce(D)*sqrt(k_)/delta(), k_)
00153     );
00154 
00155     kEqn.relax();
00156     kEqn.solve();
00157 
00158     bound(k_, k0());
00159 
00160     updateSubGridScaleFields(D);
00161 }
00162 
00163 
00164 bool dynOneEqEddy::read()
00165 {
00166     if (GenEddyVisc::read())
00167     {
00168         filter_.read(coeffDict());
00169 
00170         return true;
00171     }
00172     else
00173     {
00174         return false;
00175     }
00176 }
00177 
00178 
00179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00180 
00181 } // End namespace LESModels
00182 } // End namespace incompressible
00183 } // End namespace Foam
00184 
00185 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines