RAS models: Corrected reading of the lower-limit of epsilon and omega.
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / RAS / RASModel / RASModel.C
blob99e88408867b52983816168644967d9b03e8545f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 \*---------------------------------------------------------------------------*/
27 #include "RASModel.H"
28 #include "wallFvPatch.H"
29 #include "addToRunTimeSelectionTable.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
35 namespace incompressible
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(RASModel, 0);
41 defineRunTimeSelectionTable(RASModel, dictionary);
42 addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
44 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
46 void RASModel::printCoeffs()
48     if (printCoeffs_)
49     {
50         Info<< type() << "Coeffs" << coeffDict_ << endl;
51     }
55 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
57 RASModel::RASModel
59     const word& type,
60     const volVectorField& U,
61     const surfaceScalarField& phi,
62     transportModel& lamTransportModel
65     turbulenceModel(U, phi, lamTransportModel),
67     IOdictionary
68     (
69         IOobject
70         (
71             "RASProperties",
72             U.time().constant(),
73             U.db(),
74             IOobject::MUST_READ,
75             IOobject::NO_WRITE
76         )
77     ),
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),
89     y_(mesh_)
91     // Force the construction of the mesh deltaCoeffs which may be needed
92     // for the construction of the derived models and BCs
93     mesh_.deltaCoeffs();
97 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
99 autoPtr<RASModel> RASModel::New
101     const volVectorField& U,
102     const surfaceScalarField& phi,
103     transportModel& transport
106     word modelName;
108     // Enclose the creation of the dictionary to ensure it is deleted
109     // before the turbulenceModel is created otherwise the dictionary is
110     // entered in the database twice
111     {
112         IOdictionary dict
113         (
114             IOobject
115             (
116                 "RASProperties",
117                 U.time().constant(),
118                 U.db(),
119                 IOobject::MUST_READ,
120                 IOobject::NO_WRITE
121             )
122         );
124         dict.lookup("RASModel") >> modelName;
125     }
127     Info<< "Selecting RAS turbulence model " << modelName << endl;
129     dictionaryConstructorTable::iterator cstrIter =
130         dictionaryConstructorTablePtr_->find(modelName);
132     if (cstrIter == dictionaryConstructorTablePtr_->end())
133     {
134         FatalErrorIn
135         (
136             "RASModel::New(const volVectorField&, "
137             "const surfaceScalarField&, transportModel&)"
138         )   << "Unknown RASModel type " << modelName
139             << endl << endl
140             << "Valid RASModel types are :" << endl
141             << dictionaryConstructorTablePtr_->toc()
142             << exit(FatalError);
143     }
145     return autoPtr<RASModel>(cstrIter()(U, phi, transport));
149 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
152 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
154     scalar ypl = 11.0;
156     for (int i=0; i<10; i++)
157     {
158         ypl = log(max(E*ypl, 1))/kappa;
159     }
161     return ypl;
165 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
167     const fvPatch& curPatch = mesh_.boundary()[patchNo];
169     tmp<scalarField> tYp(new scalarField(curPatch.size()));
170     scalarField& Yp = tYp();
172     if (isA<wallFvPatch>(curPatch))
173     {
174         Yp = pow(Cmu, 0.25)
175             *y_[patchNo]
176             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
177             /nu().boundaryField()[patchNo];
178     }
179     else
180     {
181         WarningIn
182         (
183             "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
184         )   << "Patch " << patchNo << " is not a wall. Returning null field"
185             << nl << endl;
187         Yp.setSize(0);
188     }
190     return tYp;
194 void RASModel::correct()
196     turbulenceModel::correct();
198     if (turbulence_ && mesh_.changing())
199     {
200         y_.correct();
201     }
205 bool RASModel::read()
207     if (regIOobject::read())
208     {
209         lookup("turbulence") >> turbulence_;
211         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
212         {
213             coeffDict_ <<= *dictPtr;
214         }
216         k0_.readIfPresent(*this);
217         epsilon0_.readIfPresent(*this);
218         epsilonSmall_.readIfPresent(*this);
219         omega0_.readIfPresent(*this);
220         omegaSmall_.readIfPresent(*this);
222         return true;
223     }
224     else
225     {
226         return false;
227     }
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 } // End namespace incompressible
234 } // End namespace Foam
236 // ************************************************************************* //