ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / preProcessing / applyWallFunctionBoundaryConditions / applyWallFunctionBoundaryConditions.C
blob9b89af1c4de67beb66dfebe184fdf66982234f9b
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 Application
25     applyWallFunctionBounaryConditions
27 Description
28     Updates OpenFOAM RAS cases to use the new (v1.6) wall function framework
30     Attempts to determine whether case is compressible or incompressible, or
31     can be supplied with -compressible command line argument
33 \*---------------------------------------------------------------------------*/
35 #include "argList.H"
36 #include "fvMesh.H"
37 #include "Time.H"
38 #include "volFields.H"
39 #include "surfaceFields.H"
41 #include "wallPolyPatch.H"
43 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
44 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
45 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions/nutkWallFunction/nutkWallFunctionFvPatchScalarField.H"
46 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
48 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
49 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
50 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions/mutkWallFunction/mutkWallFunctionFvPatchScalarField.H"
51 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
53 using namespace Foam;
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 bool caseIsCompressible(const fvMesh& mesh)
59     // Attempt flux field
60     IOobject phiHeader
61     (
62         "phi",
63         mesh.time().timeName(),
64         mesh,
65         IOobject::MUST_READ,
66         IOobject::NO_WRITE
67     );
69     if (phiHeader.headerOk())
70     {
71         surfaceScalarField phi(phiHeader, mesh);
72         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
73         {
74             return true;
75         }
76     }
78     // Attempt density field
79     IOobject rhoHeader
80     (
81         "rho",
82         mesh.time().timeName(),
83         mesh,
84         IOobject::MUST_READ,
85         IOobject::NO_WRITE
86     );
88     if (rhoHeader.headerOk())
89     {
90         volScalarField rho(rhoHeader, mesh);
91         if (rho.dimensions() == dimDensity)
92         {
93             return true;
94         }
95     }
97     // Attempt pressure field
98     IOobject pHeader
99     (
100         "p",
101         mesh.time().timeName(),
102         mesh,
103         IOobject::MUST_READ,
104         IOobject::NO_WRITE
105     );
107     if (pHeader.headerOk())
108     {
109         volScalarField p(pHeader, mesh);
110         if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
111         {
112             return true;
113         }
114     }
116     // If none of the above are true, assume that the case is incompressible
117     return false;
121 void createVolScalarField
123     const fvMesh& mesh,
124     const word& fieldName,
125     const dimensionSet& dims
128     IOobject fieldHeader
129     (
130         fieldName,
131         mesh.time().timeName(),
132         mesh,
133         IOobject::MUST_READ,
134         IOobject::NO_WRITE
135     );
137     if (!fieldHeader.headerOk())
138     {
139         Info<< "Creating field " << fieldName << nl << endl;
141         volScalarField field
142         (
143             IOobject
144             (
145                 fieldName,
146                 mesh.time().timeName(),
147                 mesh,
148                 IOobject::NO_READ,
149                 IOobject::NO_WRITE
150             ),
151             mesh,
152             dimensionedScalar("zero", dims, 0.0)
153         );
155         field.write();
156     }
160 void replaceBoundaryType
162     const fvMesh& mesh,
163     const word& fieldName,
164     const word& boundaryType,
165     const string& boundaryValue
168     IOobject header
169     (
170         fieldName,
171         mesh.time().timeName(),
172         mesh,
173         IOobject::MUST_READ,
174         IOobject::NO_WRITE
175     );
177     if (!header.headerOk())
178     {
179         return;
180     }
182     Info<< "Updating boundary types for field " << header.name() << endl;
184     const word oldTypeName = IOdictionary::typeName;
185     const_cast<word&>(IOdictionary::typeName) = word::null;
187     IOdictionary dict(header);
189     const_cast<word&>(IOdictionary::typeName) = oldTypeName;
190     const_cast<word&>(dict.type()) = dict.headerClassName();
192     // Make a backup of the old file
193     if (mvBak(dict.objectPath(), "old"))
194     {
195         Info<< "    Backup original file to "
196             << (dict.objectPath() + ".old") << endl;
197     }
199     // Loop through boundary patches and update
200     const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
201     dictionary& boundaryDict = dict.subDict("boundaryField");
202     forAll(bMesh, patchI)
203     {
204         if (isA<wallPolyPatch>(bMesh[patchI]))
205         {
206             word patchName = bMesh[patchI].name();
207             dictionary& oldPatch = boundaryDict.subDict(patchName);
209             dictionary newPatch(dictionary::null);
210             newPatch.add("type", boundaryType);
211             newPatch.add("value", ("uniform " + boundaryValue).c_str());
213             oldPatch = newPatch;
214         }
215     }
217     Info<< "    writing updated " << dict.name() << nl << endl;
218     dict.regIOobject::write();
222 void updateCompressibleCase(const fvMesh& mesh)
224     Info<< "Case treated as compressible" << nl << endl;
225     createVolScalarField
226     (
227         mesh,
228         "mut",
229         dimArea/dimTime*dimDensity
230     );
231     replaceBoundaryType
232     (
233         mesh,
234         "mut",
235         compressible::RASModels::mutkWallFunctionFvPatchScalarField::typeName,
236         "0"
237     );
238     replaceBoundaryType
239     (
240         mesh,
241         "epsilon",
242         compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
243             typeName,
244         "0"
245     );
246     replaceBoundaryType
247     (
248         mesh,
249         "omega",
250         compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
251         "0"
252     );
253     replaceBoundaryType
254     (
255         mesh,
256         "k",
257         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
258         "0"
259     );
260     replaceBoundaryType
261     (
262         mesh,
263         "q",
264         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
265         "0"
266     );
267     replaceBoundaryType
268     (
269         mesh,
270         "R",
271         compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
272             typeName,
273         "(0 0 0 0 0 0)"
274     );
278 void updateIncompressibleCase(const fvMesh& mesh)
280     Info<< "Case treated as incompressible" << nl << endl;
281     createVolScalarField(mesh, "nut", dimArea/dimTime);
283     replaceBoundaryType
284     (
285         mesh,
286         "nut",
287         incompressible::RASModels::nutkWallFunctionFvPatchScalarField::typeName,
288         "0"
289     );
290     replaceBoundaryType
291     (
292         mesh,
293         "epsilon",
294         incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
295             typeName,
296         "0"
297     );
298     replaceBoundaryType
299     (
300         mesh,
301         "omega",
302         incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
303             typeName,
304         "0"
305     );
306     replaceBoundaryType
307     (
308         mesh,
309         "k",
310         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
311             typeName,
312         "0"
313     );
314     replaceBoundaryType
315     (
316         mesh,
317         "q",
318         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
319             typeName,
320         "0"
321     );
322     replaceBoundaryType
323     (
324         mesh,
325         "R",
326         incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
327             typeName,
328         "(0 0 0 0 0 0)"
329     );
333 int main(int argc, char *argv[])
335     #include "addTimeOptions.H"
336     argList::addBoolOption
337     (
338         "compressible",
339         "force use of compressible wall functions. default is auto-detect."
340     );
342     #include "setRootCase.H"
343     #include "createTime.H"
344     #include "createMesh.H"
346     const bool compressible = args.optionFound("compressible");
348     Info<< "Updating turbulence fields to operate using new run time "
349         << "selectable" << nl << "wall functions"
350         << nl << endl;
352     if (compressible || caseIsCompressible(mesh))
353     {
354         updateCompressibleCase(mesh);
355     }
356     else
357     {
358         updateIncompressibleCase(mesh);
359     }
361     Info<< "End\n" << endl;
363     return 0;
367 // ************************************************************************* //