ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / LESModel / LESModel.C
blob077518704c46c7d66f90844ae2d318e68fd00d55
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 "LESModel.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace incompressible
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(LESModel, 0);
39 defineRunTimeSelectionTable(LESModel, dictionary);
40 addToRunTimeSelectionTable(turbulenceModel, LESModel, turbulenceModel);
42 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
44 void LESModel::printCoeffs()
46     if (printCoeffs_)
47     {
48         Info<< type() << "Coeffs" << coeffDict_ << endl;
49     }
53 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
55 LESModel::LESModel
57     const word& type,
58     const volVectorField& U,
59     const surfaceScalarField& phi,
60     transportModel& transport,
61     const word& turbulenceModelName
64     turbulenceModel(U, phi, transport, turbulenceModelName),
66     IOdictionary
67     (
68         IOobject
69         (
70             "LESProperties",
71             U.time().constant(),
72             U.db(),
73             IOobject::MUST_READ_IF_MODIFIED,
74             IOobject::NO_WRITE
75         )
76     ),
78     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
79     coeffDict_(subOrEmptyDict(type + "Coeffs")),
81     kMin_("kMin", sqr(dimVelocity), SMALL),
82     delta_(LESdelta::New("delta", U.mesh(), *this))
84     kMin_.readIfPresent(*this);
86     // Force the construction of the mesh deltaCoeffs which may be needed
87     // for the construction of the derived models and BCs
88     mesh_.deltaCoeffs();
92 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
94 autoPtr<LESModel> LESModel::New
96     const volVectorField& U,
97     const surfaceScalarField& phi,
98     transportModel& transport,
99     const word& turbulenceModelName
102     // get model name, but do not register the dictionary
103     // otherwise it is registered in the database twice
104     const word modelType
105     (
106         IOdictionary
107         (
108             IOobject
109             (
110                 "LESProperties",
111                 U.time().constant(),
112                 U.db(),
113                 IOobject::MUST_READ_IF_MODIFIED,
114                 IOobject::NO_WRITE,
115                 false
116             )
117         ).lookup("LESModel")
118     );
120     Info<< "Selecting LES turbulence model " << modelType << endl;
122     dictionaryConstructorTable::iterator cstrIter =
123         dictionaryConstructorTablePtr_->find(modelType);
125     if (cstrIter == dictionaryConstructorTablePtr_->end())
126     {
127         FatalErrorIn
128         (
129             "LESModel::New"
130             "("
131                 "const volVectorField&, "
132                 "const surfaceScalarField& ,"
133                 "transportModel&, "
134                 "const word&"
135             ")"
136         )   << "Unknown LESModel type "
137             << modelType << nl << nl
138             << "Valid LESModel types:" << endl
139             << dictionaryConstructorTablePtr_->sortedToc()
140             << exit(FatalError);
141     }
143     return autoPtr<LESModel>
144     (
145         cstrIter()(U, phi, transport, turbulenceModelName)
146     );
150 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
152 void LESModel::correct(const tmp<volTensorField>&)
154     turbulenceModel::correct();
155     delta_().correct();
159 void LESModel::correct()
161     correct(fvc::grad(U_));
165 bool LESModel::read()
167     //if (regIOobject::read())
169     // Bit of trickery : we are both IOdictionary ('RASProperties') and
170     // an regIOobject from the turbulenceModel level. Problem is to distinguish
171     // between the two - we only want to reread the IOdictionary.
172     
173     bool ok = IOdictionary::readData
174     (
175         IOdictionary::readStream
176         (
177             IOdictionary::type()
178         )
179     );
180     IOdictionary::close();
182     if (ok)
183     {
184         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
185         {
186             coeffDict_ <<= *dictPtr;
187         }
189         delta_().read(*this);
191         kMin_.readIfPresent(*this);
193         return true;
194     }
195     else
196     {
197         return false;
198     }
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 } // End namespace incompressible
205 } // End namespace Foam
207 // ************************************************************************* //