ENH: RASModel.C: clipping input to log
[OpenFOAM-1.7.x.git] / src / turbulenceModels / compressible / RAS / RASModel / RASModel.C
blobad35e32a96c0fc6aaf1cf109ea5db79de1d22f05
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
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
19     for more details.
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 \*---------------------------------------------------------------------------*/
26 #include "RASModel.H"
27 #include "wallFvPatch.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     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
99     mesh_.deltaCoeffs();
103 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
105 autoPtr<RASModel> RASModel::New
107     const volScalarField& rho,
108     const volVectorField& U,
109     const surfaceScalarField& phi,
110     const basicThermo& thermophysicalModel
113     word modelName;
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
118     {
119         IOdictionary dict
120         (
121             IOobject
122             (
123                 "RASProperties",
124                 U.time().constant(),
125                 U.db(),
126                 IOobject::MUST_READ,
127                 IOobject::NO_WRITE
128             )
129         );
131         dict.lookup("RASModel") >> modelName;
132     }
134     Info<< "Selecting RAS turbulence model " << modelName << endl;
136     dictionaryConstructorTable::iterator cstrIter =
137         dictionaryConstructorTablePtr_->find(modelName);
139     if (cstrIter == dictionaryConstructorTablePtr_->end())
140     {
141         FatalErrorIn
142         (
143             "RASModel::New(const volScalarField&, "
144             "const volVectorField&, const surfaceScalarField&, "
145             "basicThermo&)"
146         )   << "Unknown RASModel type " << modelName
147             << endl << endl
148             << "Valid RASModel types are :" << endl
149             << dictionaryConstructorTablePtr_->sortedToc()
150             << exit(FatalError);
151     }
153     return autoPtr<RASModel>
154     (
155         cstrIter()(rho, U, phi, thermophysicalModel)
156     );
160 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
162 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
164     scalar ypl = 11.0;
166     for (int i=0; i<10; i++)
167     {
168         ypl = log(max(E*ypl, 1))/kappa;
169     }
171     return ypl;
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))
183     {
184         Yp = pow(Cmu, 0.25)
185             *y_[patchNo]
186             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
187            /(
188                 mu().boundaryField()[patchNo].patchInternalField()
189                /rho_.boundaryField()[patchNo]
190             );
191     }
192     else
193     {
194         WarningIn
195         (
196             "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
197         )   << "Patch " << patchNo << " is not a wall. Returning null field"
198             << nl << endl;
200         Yp.setSize(0);
201     }
203     return tYp;
207 void RASModel::correct()
209     if (mesh_.changing())
210     {
211         y_.correct();
212     }
216 bool RASModel::read()
218     if (regIOobject::read())
219     {
220         lookup("turbulence") >> turbulence_;
222         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
223         {
224             coeffDict_ <<= *dictPtr;
225         }
227         k0_.readIfPresent(*this);
228         epsilon0_.readIfPresent(*this);
229         epsilonSmall_.readIfPresent(*this);
230         omega0_.readIfPresent(*this);
231         omegaSmall_.readIfPresent(*this);
233         return true;
234     }
235     else
236     {
237         return false;
238     }
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 } // End namespace compressible
245 } // End namespace Foam
247 // ************************************************************************* //