1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(RASModel, 0);
39 defineRunTimeSelectionTable(RASModel, dictionary);
40 addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
42 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 void RASModel::printCoeffs()
48 Info<< type() << "Coeffs" << coeffDict_ << endl;
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 const volVectorField& U,
59 const surfaceScalarField& phi,
60 transportModel& lamTransportModel
63 turbulenceModel(U, phi, lamTransportModel),
77 turbulence_(lookup("turbulence")),
78 printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
79 coeffDict_(subOrEmptyDict(type + "Coeffs")),
81 k0_("k0", dimVelocity*dimVelocity, SMALL),
82 epsilon0_("epsilon0", k0_.dimensions()/dimTime, SMALL),
83 epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
84 omega0_("omega0", dimless/dimTime, SMALL),
85 omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
86 nuRatio_(lookupOrDefault<scalar>("nuRatio", 1e4)),
90 // Force the construction of the mesh deltaCoeffs which may be needed
91 // for the construction of the derived models and BCs
96 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
98 autoPtr<RASModel> RASModel::New
100 const volVectorField& U,
101 const surfaceScalarField& phi,
102 transportModel& transport
107 // Enclose the creation of the dictionary to ensure it is deleted
108 // before the turbulenceModel is created otherwise the dictionary is
109 // entered in the database twice
123 dict.lookup("RASModel") >> modelName;
126 Info<< "Selecting RAS turbulence model " << modelName << endl;
128 dictionaryConstructorTable::iterator cstrIter =
129 dictionaryConstructorTablePtr_->find(modelName);
131 if (cstrIter == dictionaryConstructorTablePtr_->end())
135 "RASModel::New(const volVectorField&, "
136 "const surfaceScalarField&, transportModel&)"
137 ) << "Unknown RASModel type " << modelName
139 << "Valid RASModel types are :" << endl
140 << dictionaryConstructorTablePtr_->sortedToc()
144 return autoPtr<RASModel>(cstrIter()(U, phi, transport));
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
155 for (int i=0; i<10; i++)
157 ypl = log(max(E*ypl, 1))/kappa;
164 tmp<volScalarField> RASModel::nuEff() const
166 tmp<volScalarField> tnuEff
168 new volScalarField("nuEff", nut() + nu())
172 tnuEff().internalField() =
173 Foam::min(tnuEff().internalField(), nuRatio_*nu().internalField());
179 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
181 const fvPatch& curPatch = mesh_.boundary()[patchNo];
183 tmp<scalarField> tYp(new scalarField(curPatch.size()));
184 scalarField& Yp = tYp();
186 if (curPatch.isWall())
190 *sqrt(k()().boundaryField()[patchNo].patchInternalField())
191 /nu().boundaryField()[patchNo];
197 "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
198 ) << "Patch " << patchNo << " is not a wall. Returning null field"
208 void RASModel::correct()
210 turbulenceModel::correct();
212 if (turbulence_ && mesh_.changing())
219 bool RASModel::read()
221 if (regIOobject::read())
223 lookup("turbulence") >> turbulence_;
225 if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
227 coeffDict_ <<= *dictPtr;
230 k0_.readIfPresent(*this);
231 epsilon0_.readIfPresent(*this);
232 epsilonSmall_.readIfPresent(*this);
233 omega0_.readIfPresent(*this);
234 omegaSmall_.readIfPresent(*this);
235 readIfPresent("nuRatio", nuRatio_);
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 } // End namespace incompressible
249 } // End namespace Foam
251 // ************************************************************************* //