RAS models: Corrected reading of the lower-limit of epsilon and omega.
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / RAS / RASModel / RASModel.C
blobd20dce50251856fba314adbb55525ce046dc2bdb
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 compressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
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 volScalarField& rho,
61     const volVectorField& U,
62     const surfaceScalarField& phi,
63     const basicThermo& thermophysicalModel
66     turbulenceModel(rho, U, phi, thermophysicalModel),
68     IOdictionary
69     (
70         IOobject
71         (
72             "RASProperties",
73             U.time().constant(),
74             U.db(),
75             IOobject::MUST_READ,
76             IOobject::NO_WRITE
77         )
78     ),
80     turbulence_(lookup("turbulence")),
81     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
82     coeffDict_(subOrEmptyDict(type + "Coeffs")),
84     k0_("k0", dimVelocity*dimVelocity, SMALL),
85     epsilon0_("epsilon0", k0_.dimensions()/dimTime, SMALL),
86     epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
87     omega0_("omega0", dimless/dimTime, SMALL),
88     omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
90     y_(mesh_)
92     // Force the construction of the mesh deltaCoeffs which may be needed
93     // for the construction of the derived models and BCs
94     mesh_.deltaCoeffs();
98 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
100 autoPtr<RASModel> RASModel::New
102     const volScalarField& rho,
103     const volVectorField& U,
104     const surfaceScalarField& phi,
105     const basicThermo& thermophysicalModel
108     word modelName;
110     // Enclose the creation of the dictionary to ensure it is deleted
111     // before the turbulenceModel is created otherwise the dictionary is
112     // entered in the database twice
113     {
114         IOdictionary dict
115         (
116             IOobject
117             (
118                 "RASProperties",
119                 U.time().constant(),
120                 U.db(),
121                 IOobject::MUST_READ,
122                 IOobject::NO_WRITE
123             )
124         );
126         dict.lookup("RASModel") >> modelName;
127     }
129     Info<< "Selecting RAS turbulence model " << modelName << endl;
131     dictionaryConstructorTable::iterator cstrIter =
132         dictionaryConstructorTablePtr_->find(modelName);
134     if (cstrIter == dictionaryConstructorTablePtr_->end())
135     {
136         FatalErrorIn
137         (
138             "RASModel::New(const volScalarField&, "
139             "const volVectorField&, const surfaceScalarField&, "
140             "basicThermo&)"
141         )   << "Unknown RASModel type " << modelName
142             << endl << endl
143             << "Valid RASModel types are :" << endl
144             << dictionaryConstructorTablePtr_->toc()
145             << exit(FatalError);
146     }
148     return autoPtr<RASModel>
149     (
150         cstrIter()(rho, U, phi, thermophysicalModel)
151     );
155 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
157 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
159     scalar ypl = 11.0;
161     for (int i=0; i<10; i++)
162     {
163         ypl = log(E*ypl)/kappa;
164     }
166     return ypl;
170 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
172     const fvPatch& curPatch = mesh_.boundary()[patchNo];
174     tmp<scalarField> tYp(new scalarField(curPatch.size()));
175     scalarField& Yp = tYp();
177     if (isA<wallFvPatch>(curPatch))
178     {
179         Yp = pow(Cmu, 0.25)
180             *y_[patchNo]
181             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
182            /(
183                 mu().boundaryField()[patchNo].patchInternalField()
184                /rho_.boundaryField()[patchNo]
185             );
186     }
187     else
188     {
189         WarningIn
190         (
191             "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
192         )   << "Patch " << patchNo << " is not a wall. Returning null field"
193             << nl << endl;
195         Yp.setSize(0);
196     }
198     return tYp;
202 void RASModel::correct()
204     if (mesh_.changing())
205     {
206         y_.correct();
207     }
211 bool RASModel::read()
213     if (regIOobject::read())
214     {
215         lookup("turbulence") >> turbulence_;
217         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
218         {
219             coeffDict_ <<= *dictPtr;
220         }
222         k0_.readIfPresent(*this);
223         epsilon0_.readIfPresent(*this);
224         epsilonSmall_.readIfPresent(*this);
225         omega0_.readIfPresent(*this);
226         omegaSmall_.readIfPresent(*this);
228         return true;
229     }
230     else
231     {
232         return false;
233     }
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 } // End namespace compressible
240 } // End namespace Foam
242 // ************************************************************************* //