1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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()
48 Info<< type() << "Coeffs" << coeffDict_ << endl;
53 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
58 const volScalarField& rho,
59 const volVectorField& U,
60 const surfaceScalarField& phi,
61 const basicThermo& thermoPhysicalModel,
62 const word& turbulenceModelName
65 turbulenceModel(rho, U, phi, thermoPhysicalModel, turbulenceModelName),
74 IOobject::MUST_READ_IF_MODIFIED,
79 printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
80 coeffDict_(subOrEmptyDict(type + "Coeffs")),
82 kMin_("kMin", sqr(dimVelocity), SMALL),
84 delta_(LESdelta::New("delta", U.mesh(), *this))
86 kMin_.readIfPresent(*this);
88 // Force the construction of the mesh deltaCoeffs which may be needed
89 // for the construction of the derived models and BCs
94 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
96 autoPtr<LESModel> LESModel::New
98 const volScalarField& rho,
99 const volVectorField& U,
100 const surfaceScalarField& phi,
101 const basicThermo& thermoPhysicalModel,
102 const word& turbulenceModelName
105 // get model name, but do not register the dictionary
106 // otherwise it is registered in the database twice
116 IOobject::MUST_READ_IF_MODIFIED,
123 Info<< "Selecting LES turbulence model " << modelType << endl;
125 dictionaryConstructorTable::iterator cstrIter =
126 dictionaryConstructorTablePtr_->find(modelType);
128 if (cstrIter == dictionaryConstructorTablePtr_->end())
134 "const volScalarField&, "
135 "const volVectorField&, "
136 "const surfaceScalarField&, "
137 "const basicThermo&, "
140 ) << "Unknown LESModel type "
141 << modelType << nl << nl
142 << "Valid LESModel types:" << endl
143 << dictionaryConstructorTablePtr_->sortedToc()
147 return autoPtr<LESModel>
149 cstrIter()(rho, U, phi, thermoPhysicalModel, turbulenceModelName)
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 void LESModel::correct(const tmp<volTensorField>&)
162 void LESModel::correct()
164 correct(fvc::grad(U_));
168 bool LESModel::read()
170 // Bit of trickery : we are both IOdictionary ('RASProperties') and
171 // an regIOobject (from the turbulenceModel). Problem is to distinguish
172 // between the two - we only want to reread the IOdictionary.
174 bool ok = IOdictionary::readData
176 IOdictionary::readStream
181 IOdictionary::close();
185 if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
187 coeffDict_ <<= *dictPtr;
190 kMin_.readIfPresent(*this);
192 delta_().read(*this);
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 } // End namespace compressible
206 } // End namespace Foam
208 // ************************************************************************* //