BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / RASModel / RASModel.C
blob99697c346097eefa8a630bbf0365b207a888f586
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 incompressible
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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 volVectorField& U,
60     const surfaceScalarField& phi,
61     transportModel& transport,
62     const word& turbulenceModelName
65     turbulenceModel(U, phi, transport, turbulenceModelName),
67     IOdictionary
68     (
69         IOobject
70         (
71             "RASProperties",
72             U.time().constant(),
73             U.db(),
74             IOobject::MUST_READ_IF_MODIFIED,
75             IOobject::NO_WRITE
76         )
77     ),
79     turbulence_(lookup("turbulence")),
80     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
81     coeffDict_(subOrEmptyDict(type + "Coeffs")),
83     kMin_("kMin", sqr(dimVelocity), SMALL),
84     epsilonMin_("epsilonMin", kMin_.dimensions()/dimTime, SMALL),
85     omegaMin_("omegaMin", dimless/dimTime, SMALL),
87     y_(mesh_)
89     kMin_.readIfPresent(*this);
90     epsilonMin_.readIfPresent(*this);
91     omegaMin_.readIfPresent(*this);
93     // Force the construction of the mesh deltaCoeffs which may be needed
94     // for the construction of the derived models and BCs
95     mesh_.deltaCoeffs();
99 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
101 autoPtr<RASModel> RASModel::New
103     const volVectorField& U,
104     const surfaceScalarField& phi,
105     transportModel& transport,
106     const word& turbulenceModelName
109     // get model name, but do not register the dictionary
110     // otherwise it is registered in the database twice
111     const word modelType
112     (
113         IOdictionary
114         (
115             IOobject
116             (
117                 "RASProperties",
118                 U.time().constant(),
119                 U.db(),
120                 IOobject::MUST_READ_IF_MODIFIED,
121                 IOobject::NO_WRITE,
122                 false
123             )
124         ).lookup("RASModel")
125     );
127     Info<< "Selecting RAS turbulence model " << modelType << endl;
129     dictionaryConstructorTable::iterator cstrIter =
130         dictionaryConstructorTablePtr_->find(modelType);
132     if (cstrIter == dictionaryConstructorTablePtr_->end())
133     {
134         FatalErrorIn
135         (
136             "RASModel::New"
137             "("
138                 "const volVectorField&, "
139                 "const surfaceScalarField&, "
140                 "transportModel&, "
141                 "const word&"
142             ")"
143         )   << "Unknown RASModel type "
144             << modelType << nl << nl
145             << "Valid RASModel types:" << endl
146             << dictionaryConstructorTablePtr_->sortedToc()
147             << exit(FatalError);
148     }
150     return autoPtr<RASModel>
151     (
152         cstrIter()(U, phi, transport, turbulenceModelName)
153     );
157 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
160 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
162     scalar ypl = 11.0;
164     for (int i=0; i<10; i++)
165     {
166         ypl = log(max(E*ypl, 1))/kappa;
167     }
169     return ypl;
173 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
175     const fvPatch& curPatch = mesh_.boundary()[patchNo];
177     tmp<scalarField> tYp(new scalarField(curPatch.size()));
178     scalarField& Yp = tYp();
180     if (isA<wallFvPatch>(curPatch))
181     {
182         Yp = pow025(Cmu)
183             *y_[patchNo]
184             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
185             /nu()().boundaryField()[patchNo];
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     turbulenceModel::correct();
206     if (turbulence_ && mesh_.changing())
207     {
208         y_.correct();
209     }
213 bool RASModel::read()
215     //if (regIOobject::read())
217     // Bit of trickery : we are both IOdictionary ('RASProperties') and
218     // an regIOobject from the turbulenceModel level. Problem is to distinguish
219     // between the two - we only want to reread the IOdictionary.
220     
221     bool ok = IOdictionary::readData
222     (
223         IOdictionary::readStream
224         (
225             IOdictionary::type()
226         )
227     );
228     IOdictionary::close();
230     if (ok)
231     {
232         lookup("turbulence") >> turbulence_;
234         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
235         {
236             coeffDict_ <<= *dictPtr;
237         }
239         kMin_.readIfPresent(*this);
240         epsilonMin_.readIfPresent(*this);
241         omegaMin_.readIfPresent(*this);
243         return true;
244     }
245     else
246     {
247         return false;
248     }
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 } // End namespace incompressible
255 } // End namespace Foam
257 // ************************************************************************* //