initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / LES / LESModel / LESModel.C
blob5602e44d7f8ea3b6e3e8746b53085b2d7faa2315
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 "LESModel.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace compressible
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(LESModel, 0);
40 defineRunTimeSelectionTable(LESModel, dictionary);
41 addToRunTimeSelectionTable(turbulenceModel, LESModel, turbulenceModel);
43 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
45 void LESModel::printCoeffs()
47     if (printCoeffs_)
48     {
49         Info<< type() << "Coeffs" << coeffDict_ << endl;
50     }
54 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
56 LESModel::LESModel
58     const word& type,
59     const volScalarField& rho,
60     const volVectorField& U,
61     const surfaceScalarField& phi,
62     const basicThermo& thermoPhysicalModel
65     turbulenceModel(rho, U, phi, thermoPhysicalModel),
67     IOdictionary
68     (
69         IOobject
70         (
71             "LESProperties",
72             U.time().constant(),
73             U.db(),
74             IOobject::MUST_READ,
75             IOobject::NO_WRITE
76         )
77     ),
79     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
80     coeffDict_(subDictPtr(type + "Coeffs")),
82     k0_("k0", dimVelocity*dimVelocity, SMALL),
84     delta_(LESdelta::New("delta", U.mesh(), *this))
86     readIfPresent("k0", k0_);
90 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
92 autoPtr<LESModel> LESModel::New
94     const volScalarField& rho,
95     const volVectorField& U,
96     const surfaceScalarField& phi,
97     const basicThermo& thermoPhysicalModel
100     word modelName;
102     // Enclose the creation of the dictionary to ensure it is deleted
103     // before the turbulenceModel is created otherwise the dictionary is
104     // entered in the database twice
105     {
106         IOdictionary dict
107         (
108             IOobject
109             (
110                 "LESProperties",
111                 U.time().constant(),
112                 U.db(),
113                 IOobject::MUST_READ,
114                 IOobject::NO_WRITE
115             )
116         );
118         dict.lookup("LESModel") >> modelName;
119     }
121     Info<< "Selecting LES turbulence model " << modelName << endl;
123     dictionaryConstructorTable::iterator cstrIter =
124         dictionaryConstructorTablePtr_->find(modelName);
126     if (cstrIter == dictionaryConstructorTablePtr_->end())
127     {
128         FatalErrorIn
129         (
130             "LESModel::New(const volVectorField& U, const "
131             "surfaceScalarField& phi, const basicThermo&)"
132         )   << "Unknown LESModel type " << modelName
133             << endl << endl
134             << "Valid LESModel types are :" << endl
135             << dictionaryConstructorTablePtr_->toc()
136             << exit(FatalError);
137     }
139     return autoPtr<LESModel>(cstrIter()(rho, U, phi, thermoPhysicalModel));
143 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
145 void LESModel::correct(const tmp<volTensorField>&)
147     delta_().correct();
151 void LESModel::correct()
153     correct(fvc::grad(U_));
157 bool LESModel::read()
159     if (regIOobject::read())
160     {
161         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
162         {
163             coeffDict_ <<= *dictPtr;
164         }
166         readIfPresent("k0", k0_);
168         delta_().read(*this);
170         return true;
171     }
172     else
173     {
174         return false;
175     }
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 } // End namespace compressible
182 } // End namespace Foam
184 // ************************************************************************* //