Merge commit 'd3b269b7c6ffa0cd68845adfecdfb849316dba71' into nextRelease
[foam-extend-3.2.git] / src / turbulenceModels / compressible / LES / LESModel / LESModel.C
blob410e2b6a428c11224e64b2932be32f03998e7766
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     |
5     \\  /    A nd           | For copyright notice see file Copyright
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "LESModel.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace compressible
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(LESModel, 0);
39 defineRunTimeSelectionTable(LESModel, dictionary);
40 addToRunTimeSelectionTable(turbulenceModel, LESModel, turbulenceModel);
42 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
44 void LESModel::printCoeffs()
46     if (printCoeffs_)
47     {
48         Info<< type() << "Coeffs" << coeffDict_ << endl;
49     }
53 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
55 LESModel::LESModel
57     const word& type,
58     const volScalarField& rho,
59     const volVectorField& U,
60     const surfaceScalarField& phi,
61     const basicThermo& thermoPhysicalModel
64     turbulenceModel(rho, U, phi, thermoPhysicalModel),
66     IOdictionary
67     (
68         IOobject
69         (
70             "LESProperties",
71             U.time().constant(),
72             U.db(),
73             IOobject::MUST_READ,
74             IOobject::NO_WRITE
75         )
76     ),
78     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
79     coeffDict_(subOrEmptyDict(type + "Coeffs")),
81     k0_("k0", dimVelocity*dimVelocity, SMALL),
83     delta_(LESdelta::New("delta", U.mesh(), *this))
85     readIfPresent("k0", k0_);
87     // Force the construction of the mesh deltaCoeffs which may be needed
88     // for the construction of the derived models and BCs
89     mesh_.deltaCoeffs();
93 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
95 autoPtr<LESModel> LESModel::New
97     const volScalarField& rho,
98     const volVectorField& U,
99     const surfaceScalarField& phi,
100     const basicThermo& thermoPhysicalModel
103     word modelName;
105     // Enclose the creation of the dictionary to ensure it is deleted
106     // before the turbulenceModel is created otherwise the dictionary is
107     // entered in the database twice
108     {
109         IOdictionary dict
110         (
111             IOobject
112             (
113                 "LESProperties",
114                 U.time().constant(),
115                 U.db(),
116                 IOobject::MUST_READ,
117                 IOobject::NO_WRITE
118             )
119         );
121         dict.lookup("LESModel") >> modelName;
122     }
124     Info<< "Selecting LES turbulence model " << modelName << endl;
126     dictionaryConstructorTable::iterator cstrIter =
127         dictionaryConstructorTablePtr_->find(modelName);
129     if (cstrIter == dictionaryConstructorTablePtr_->end())
130     {
131         FatalErrorIn
132         (
133             "LESModel::New(const volVectorField& U, const "
134             "surfaceScalarField& phi, const basicThermo&)"
135         )   << "Unknown LESModel type " << modelName
136             << endl << endl
137             << "Valid LESModel types are :" << endl
138             << dictionaryConstructorTablePtr_->sortedToc()
139             << exit(FatalError);
140     }
142     return autoPtr<LESModel>(cstrIter()(rho, U, phi, thermoPhysicalModel));
146 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
148 void LESModel::correct(const tmp<volTensorField>&)
150     delta_().correct();
154 void LESModel::correct()
156     correct(fvc::grad(U_));
160 bool LESModel::read()
162     if (regIOobject::read())
163     {
164         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
165         {
166             coeffDict_ <<= *dictPtr;
167         }
169         readIfPresent("k0", k0_);
171         delta_().read(*this);
173         return true;
174     }
175     else
176     {
177         return false;
178     }
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 } // End namespace compressible
185 } // End namespace Foam
187 // ************************************************************************* //