changed isA<wallPolyPatch> & isA<wallFvPatch> almost everywhere
[openfoam-extend-OpenFOAM-1.6-ext.git] / src / turbulenceModels / compressible / RAS / RASModel / RASModel.C
blobd72f0c54abb357f6b6db78c79a8c62c13c9f0e99
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
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()
47     if (printCoeffs_)
48     {
49         Info<< type() << "Coeffs" << coeffDict_ << endl;
50     }
54 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
56 RASModel::RASModel
58     const word& type,
59     const volScalarField& rho,
60     const volVectorField& U,
61     const surfaceScalarField& phi,
62     const basicThermo& thermophysicalModel
65     turbulenceModel(rho, U, phi, thermophysicalModel),
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 volScalarField& rho,
102     const volVectorField& U,
103     const surfaceScalarField& phi,
104     const basicThermo& thermophysicalModel
107     word modelName;
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
112     {
113         IOdictionary dict
114         (
115             IOobject
116             (
117                 "RASProperties",
118                 U.time().constant(),
119                 U.db(),
120                 IOobject::MUST_READ,
121                 IOobject::NO_WRITE
122             )
123         );
125         dict.lookup("RASModel") >> modelName;
126     }
128     Info<< "Selecting RAS turbulence model " << modelName << endl;
130     dictionaryConstructorTable::iterator cstrIter =
131         dictionaryConstructorTablePtr_->find(modelName);
133     if (cstrIter == dictionaryConstructorTablePtr_->end())
134     {
135         FatalErrorIn
136         (
137             "RASModel::New(const volScalarField&, "
138             "const volVectorField&, const surfaceScalarField&, "
139             "basicThermo&)"
140         )   << "Unknown RASModel type " << modelName
141             << endl << endl
142             << "Valid RASModel types are :" << endl
143             << dictionaryConstructorTablePtr_->toc()
144             << exit(FatalError);
145     }
147     return autoPtr<RASModel>
148     (
149         cstrIter()(rho, U, phi, thermophysicalModel)
150     );
154 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
156 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
158     scalar ypl = 11.0;
160     for (int i=0; i<10; i++)
161     {
162         ypl = log(E*ypl)/kappa;
163     }
165     return ypl;
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())
177     {
178         Yp = pow(Cmu, 0.25)
179             *y_[patchNo]
180             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
181            /(
182                 mu().boundaryField()[patchNo].patchInternalField()
183                /rho_.boundaryField()[patchNo]
184             );
185     }
186     else
187     {
188         WarningIn
189         (
190             "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
191         )   << "Patch " << patchNo << " is not a wall. Returning null field"
192             << nl << endl;
194         Yp.setSize(0);
195     }
197     return tYp;
201 void RASModel::correct()
203     if (mesh_.changing())
204     {
205         y_.correct();
206     }
210 bool RASModel::read()
212     if (regIOobject::read())
213     {
214         lookup("turbulence") >> turbulence_;
216         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
217         {
218             coeffDict_ <<= *dictPtr;
219         }
221         k0_.readIfPresent(*this);
222         epsilon0_.readIfPresent(*this);
223         epsilonSmall_.readIfPresent(*this);
224         omega0_.readIfPresent(*this);
225         omegaSmall_.readIfPresent(*this);
227         return true;
228     }
229     else
230     {
231         return false;
232     }
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 } // End namespace compressible
239 } // End namespace Foam
241 // ************************************************************************* //