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

spectEddyVisc.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 "spectEddyVisc.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(spectEddyVisc, 0);
00041 addToRunTimeSelectionTable(LESModel, spectEddyVisc, dictionary);
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void spectEddyVisc::updateSubGridScaleFields(const volTensorField& gradU)
00046 {
00047     volScalarField Re = sqr(delta())*mag(symm(gradU))/nu();
00048     for (label i=0; i<5; i++)
00049     {
00050         nuSgs_ =
00051             nu()
00052            /(
00053                  scalar(1)
00054                - exp(-cB_*pow(nu()/(nuSgs_ + nu()), 1.0/3.0)*pow(Re, -2.0/3.0))
00055             );
00056     }
00057 
00058     nuSgs_.correctBoundaryConditions();
00059 }
00060 
00061 
00062 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00063 
00064 spectEddyVisc::spectEddyVisc
00065 (
00066     const volVectorField& U,
00067     const surfaceScalarField& phi,
00068     transportModel& transport
00069 )
00070 :
00071     LESModel(typeName, U, phi, transport),
00072     GenEddyVisc(U, phi, transport),
00073 
00074 
00075     cB_
00076     (
00077         dimensioned<scalar>::lookupOrAddToDict
00078         (
00079             "cB",
00080             coeffDict_,
00081             8.22
00082         )
00083     ),
00084     cK1_
00085     (
00086         dimensioned<scalar>::lookupOrAddToDict
00087         (
00088             "cK1",
00089             coeffDict_,
00090             0.83
00091         )
00092     ),
00093     cK2_
00094     (
00095         dimensioned<scalar>::lookupOrAddToDict
00096         (
00097             "cK2",
00098             coeffDict_,
00099             1.03
00100         )
00101     ),
00102     cK3_
00103     (
00104         dimensioned<scalar>::lookupOrAddToDict
00105         (
00106             "cK3",
00107             coeffDict_,
00108             4.75
00109         )
00110     ),
00111     cK4_
00112     (
00113         dimensioned<scalar>::lookupOrAddToDict
00114         (
00115             "cK4",
00116             coeffDict_,
00117             2.55
00118         )
00119     )
00120 {
00121     printCoeffs();
00122 
00123     updateSubGridScaleFields(fvc::grad(U));
00124 }
00125 
00126 
00127 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00128 
00129 tmp<volScalarField> spectEddyVisc::k() const
00130 {
00131     volScalarField eps = 2*nuEff()*magSqr(symm(fvc::grad(U())));
00132 
00133     return
00134         cK1_*pow(delta()*eps, 2.0/3.0)
00135        *exp(-cK2_*pow(delta(), -4.0/3.0)*nu()/pow(eps, 1.0/3.0))
00136       - cK3_*sqrt(eps*nu())
00137        *erfc(cK4_*pow(delta(), -2.0/3.0)*sqrt(nu())*pow(eps, -1.0/6.0));
00138 }
00139 
00140 
00141 void spectEddyVisc::correct(const tmp<volTensorField>& gradU)
00142 {
00143     GenEddyVisc::correct(gradU);
00144     updateSubGridScaleFields(gradU());
00145 }
00146 
00147 
00148 bool spectEddyVisc::read()
00149 {
00150     if (GenEddyVisc::read())
00151     {
00152         cB_.readIfPresent(coeffDict());
00153         cK1_.readIfPresent(coeffDict());
00154         cK2_.readIfPresent(coeffDict());
00155         cK3_.readIfPresent(coeffDict());
00156         cK4_.readIfPresent(coeffDict());
00157 
00158         return true;
00159     }
00160     else
00161     {
00162         return false;
00163     }
00164 }
00165 
00166 
00167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00168 
00169 } // End namespace LESModels
00170 } // End namespace incompressible
00171 } // End namespace Foam
00172 
00173 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines