1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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 "wallFvPatch.H"
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 k0_.readIfPresent(*this);
92 epsilon0_.readIfPresent(*this);
93 epsilonSmall_.readIfPresent(*this);
94 omega0_.readIfPresent(*this);
95 omegaSmall_.readIfPresent(*this);
97 // Force the construction of the mesh deltaCoeffs which may be needed
98 // for the construction of the derived models and BCs
103 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
105 autoPtr<RASModel> RASModel::New
107 const volScalarField& rho,
108 const volVectorField& U,
109 const surfaceScalarField& phi,
110 const basicThermo& thermophysicalModel
115 // Enclose the creation of the dictionary to ensure it is deleted
116 // before the turbulenceModel is created otherwise the dictionary is
117 // entered in the database twice
131 dict.lookup("RASModel") >> modelName;
134 Info<< "Selecting RAS turbulence model " << modelName << endl;
136 dictionaryConstructorTable::iterator cstrIter =
137 dictionaryConstructorTablePtr_->find(modelName);
139 if (cstrIter == dictionaryConstructorTablePtr_->end())
143 "RASModel::New(const volScalarField&, "
144 "const volVectorField&, const surfaceScalarField&, "
146 ) << "Unknown RASModel type " << modelName
148 << "Valid RASModel types are :" << endl
149 << dictionaryConstructorTablePtr_->sortedToc()
153 return autoPtr<RASModel>
155 cstrIter()(rho, U, phi, thermophysicalModel)
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
166 for (int i=0; i<10; i++)
168 ypl = log(max(E*ypl, 1))/kappa;
175 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
177 const fvPatch& curPatch = mesh_.boundary()[patchNo];
179 tmp<scalarField> tYp(new scalarField(curPatch.size()));
180 scalarField& Yp = tYp();
182 if (isA<wallFvPatch>(curPatch))
186 *sqrt(k()().boundaryField()[patchNo].patchInternalField())
188 mu().boundaryField()[patchNo].patchInternalField()
189 /rho_.boundaryField()[patchNo]
196 "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
197 ) << "Patch " << patchNo << " is not a wall. Returning null field"
207 void RASModel::correct()
209 if (mesh_.changing())
216 bool RASModel::read()
218 if (regIOobject::read())
220 lookup("turbulence") >> turbulence_;
222 if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
224 coeffDict_ <<= *dictPtr;
227 k0_.readIfPresent(*this);
228 epsilon0_.readIfPresent(*this);
229 epsilonSmall_.readIfPresent(*this);
230 omega0_.readIfPresent(*this);
231 omegaSmall_.readIfPresent(*this);
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 } // End namespace compressible
245 } // End namespace Foam
247 // ************************************************************************* //