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

oneEqEddy.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 "oneEqEddy.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace compressible
00034 {
00035 namespace LESModels
00036 {
00037 
00038 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(oneEqEddy, 0);
00041 addToRunTimeSelectionTable(LESModel, oneEqEddy, dictionary);
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void oneEqEddy::updateSubGridScaleFields()
00046 {
00047     muSgs_ = ck_*rho()*sqrt(k_)*delta();
00048     muSgs_.correctBoundaryConditions();
00049 
00050     alphaSgs_ = muSgs_/Prt_;
00051     alphaSgs_.correctBoundaryConditions();
00052 }
00053 
00054 
00055 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00056 
00057 oneEqEddy::oneEqEddy
00058 (
00059     const volScalarField& rho,
00060     const volVectorField& U,
00061     const surfaceScalarField& phi,
00062     const basicThermo& thermoPhysicalModel
00063 )
00064 :
00065     LESModel(typeName, rho, U, phi, thermoPhysicalModel),
00066     GenEddyVisc(rho, U, phi, thermoPhysicalModel),
00067 
00068     ck_
00069     (
00070         dimensioned<scalar>::lookupOrAddToDict
00071         (
00072             "ck",
00073             coeffDict_,
00074             0.094
00075         )
00076     )
00077 {
00078     updateSubGridScaleFields();
00079 
00080     printCoeffs();
00081 }
00082 
00083 
00084 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00085 
00086 void oneEqEddy::correct(const tmp<volTensorField>& tgradU)
00087 {
00088     const volTensorField& gradU = tgradU();
00089 
00090     GenEddyVisc::correct(gradU);
00091 
00092     volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
00093     volScalarField G = 2*muSgs_*(gradU && dev(symm(gradU)));
00094 
00095     fvScalarMatrix kEqn
00096     (
00097         fvm::ddt(rho(), k_)
00098       + fvm::div(phi(), k_)
00099       - fvm::laplacian(DkEff(), k_)
00100      ==
00101         G
00102       - fvm::SuSp(2.0/3.0*rho()*divU, k_)
00103       - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
00104     );
00105 
00106     kEqn.relax();
00107     kEqn.solve();
00108 
00109     bound(k_, k0());
00110 
00111     updateSubGridScaleFields();
00112 }
00113 
00114 
00115 bool oneEqEddy::read()
00116 {
00117     if (GenEddyVisc::read())
00118     {
00119         ck_.readIfPresent(coeffDict());
00120 
00121         return true;
00122     }
00123     else
00124     {
00125         return false;
00126     }
00127 }
00128 
00129 
00130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00131 
00132 } // End namespace LESModels
00133 } // End namespace compressible
00134 } // End namespace Foam
00135 
00136 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines