1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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 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
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 \*---------------------------------------------------------------------------*/
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace compressible
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(RASModel, 0);
40 defineRunTimeSelectionTable(RASModel, dictionary);
41 addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45 void RASModel::printCoeffs()
49 Info<< type() << "Coeffs" << coeffDict_ << endl;
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 const volScalarField& rho,
60 const volVectorField& U,
61 const surfaceScalarField& phi,
62 const basicThermo& thermophysicalModel
65 turbulenceModel(rho, U, phi, thermophysicalModel),
79 turbulence_(lookup("turbulence")),
80 printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
81 coeffDict_(subOrEmptyDict(type + "Coeffs")),
83 k0_("k0", dimVelocity*dimVelocity, SMALL),
84 epsilon0_("epsilon0", k0_.dimensions()/dimTime, SMALL),
85 epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
86 omega0_("omega0", dimless/dimTime, SMALL),
87 omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
91 // Force the construction of the mesh deltaCoeffs which may be needed
92 // for the construction of the derived models and BCs
97 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
99 autoPtr<RASModel> RASModel::New
101 const volScalarField& rho,
102 const volVectorField& U,
103 const surfaceScalarField& phi,
104 const basicThermo& thermophysicalModel
109 // Enclose the creation of the dictionary to ensure it is deleted
110 // before the turbulenceModel is created otherwise the dictionary is
111 // entered in the database twice
125 dict.lookup("RASModel") >> modelName;
128 Info<< "Selecting RAS turbulence model " << modelName << endl;
130 dictionaryConstructorTable::iterator cstrIter =
131 dictionaryConstructorTablePtr_->find(modelName);
133 if (cstrIter == dictionaryConstructorTablePtr_->end())
137 "RASModel::New(const volScalarField&, "
138 "const volVectorField&, const surfaceScalarField&, "
140 ) << "Unknown RASModel type " << modelName
142 << "Valid RASModel types are :" << endl
143 << dictionaryConstructorTablePtr_->toc()
147 return autoPtr<RASModel>
149 cstrIter()(rho, U, phi, thermophysicalModel)
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
160 for (int i=0; i<10; i++)
162 ypl = log(E*ypl)/kappa;
169 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
171 const fvPatch& curPatch = mesh_.boundary()[patchNo];
173 tmp<scalarField> tYp(new scalarField(curPatch.size()));
174 scalarField& Yp = tYp();
176 if (curPatch.isWall())
180 *sqrt(k()().boundaryField()[patchNo].patchInternalField())
182 mu().boundaryField()[patchNo].patchInternalField()
183 /rho_.boundaryField()[patchNo]
190 "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
191 ) << "Patch " << patchNo << " is not a wall. Returning null field"
201 void RASModel::correct()
203 if (mesh_.changing())
210 bool RASModel::read()
212 if (regIOobject::read())
214 lookup("turbulence") >> turbulence_;
216 if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
218 coeffDict_ <<= *dictPtr;
221 k0_.readIfPresent(*this);
222 epsilon0_.readIfPresent(*this);
223 epsilonSmall_.readIfPresent(*this);
224 omega0_.readIfPresent(*this);
225 omegaSmall_.readIfPresent(*this);
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 } // End namespace compressible
239 } // End namespace Foam
241 // ************************************************************************* //