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

LESModel.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 "LESModel.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace incompressible
00034 {
00035 
00036 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00037 
00038 defineTypeNameAndDebug(LESModel, 0);
00039 defineRunTimeSelectionTable(LESModel, dictionary);
00040 addToRunTimeSelectionTable(turbulenceModel, LESModel, turbulenceModel);
00041 
00042 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
00043 
00044 void LESModel::printCoeffs()
00045 {
00046     if (printCoeffs_)
00047     {
00048         Info<< type() << "Coeffs" << coeffDict_ << endl;
00049     }
00050 }
00051 
00052 
00053 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
00054 
00055 LESModel::LESModel
00056 (
00057     const word& type,
00058     const volVectorField& U,
00059     const surfaceScalarField& phi,
00060     transportModel& lamTransportModel
00061 )
00062 :
00063     turbulenceModel(U, phi, lamTransportModel),
00064 
00065     IOdictionary
00066     (
00067         IOobject
00068         (
00069             "LESProperties",
00070             U.time().constant(),
00071             U.db(),
00072             IOobject::MUST_READ,
00073             IOobject::NO_WRITE
00074         )
00075     ),
00076 
00077     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
00078     coeffDict_(subOrEmptyDict(type + "Coeffs")),
00079 
00080     k0_("k0", dimVelocity*dimVelocity, SMALL),
00081     delta_(LESdelta::New("delta", U.mesh(), *this))
00082 {
00083     readIfPresent("k0", k0_);
00084 
00085     // Force the construction of the mesh deltaCoeffs which may be needed
00086     // for the construction of the derived models and BCs
00087     mesh_.deltaCoeffs();
00088 }
00089 
00090 
00091 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
00092 
00093 autoPtr<LESModel> LESModel::New
00094 (
00095     const volVectorField& U,
00096     const surfaceScalarField& phi,
00097     transportModel& transport
00098 )
00099 {
00100     word modelName;
00101 
00102     // Enclose the creation of the dictionary to ensure it is deleted
00103     // before the turbulenceModel is created otherwise the dictionary is
00104     // entered in the database twice
00105     {
00106         IOdictionary dict
00107         (
00108             IOobject
00109             (
00110                 "LESProperties",
00111                 U.time().constant(),
00112                 U.db(),
00113                 IOobject::MUST_READ,
00114                 IOobject::NO_WRITE
00115             )
00116         );
00117 
00118         dict.lookup("LESModel") >> modelName;
00119     }
00120 
00121     Info<< "Selecting LES turbulence model " << modelName << endl;
00122 
00123     dictionaryConstructorTable::iterator cstrIter =
00124         dictionaryConstructorTablePtr_->find(modelName);
00125 
00126     if (cstrIter == dictionaryConstructorTablePtr_->end())
00127     {
00128         FatalErrorIn
00129         (
00130             "LESModel::New(const volVectorField& U, const "
00131             "surfaceScalarField& phi, transportModel&)"
00132         )   << "Unknown LESModel type " << modelName
00133             << endl << endl
00134             << "Valid LESModel types are :" << endl
00135             << dictionaryConstructorTablePtr_->sortedToc()
00136             << exit(FatalError);
00137     }
00138 
00139     return autoPtr<LESModel>(cstrIter()(U, phi, transport));
00140 }
00141 
00142 
00143 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00144 
00145 void LESModel::correct(const tmp<volTensorField>&)
00146 {
00147     turbulenceModel::correct();
00148     delta_().correct();
00149 }
00150 
00151 
00152 void LESModel::correct()
00153 {
00154     correct(fvc::grad(U_));
00155 }
00156 
00157 
00158 bool LESModel::read()
00159 {
00160     if (regIOobject::read())
00161     {
00162         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
00163         {
00164             coeffDict_ <<= *dictPtr;
00165         }
00166 
00167         delta_().read(*this);
00168 
00169         readIfPresent("k0", k0_);
00170 
00171         return true;
00172     }
00173     else
00174     {
00175         return false;
00176     }
00177 }
00178 
00179 
00180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00181 
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