ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / applications / utilities / preProcessing / applyWallFunctionBoundaryConditions / applyWallFunctionBoundaryConditions.C
blob32495c1cd1b28df5eef5dc2a86be0baaadaa23a0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2007 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 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/nutWallFunction/nutWallFunctionFvPatchScalarField.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/mutWallFunction/mutWallFunctionFvPatchScalarField.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 field
193     word backupName(dict.name() + ".old");
194     Info<< "    copying " << dict.name() << " to "
195         << backupName << endl;
196     IOdictionary dictOld = dict;
197     dictOld.rename(backupName);
198     dictOld.regIOobject::write();
200     // Loop through boundary patches and update
201     const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
202     dictionary& boundaryDict = dict.subDict("boundaryField");
203     forAll(bMesh, patchI)
204     {
205         if (isA<wallPolyPatch>(bMesh[patchI]))
206         {
207             word patchName = bMesh[patchI].name();
208             dictionary& oldPatch = boundaryDict.subDict(patchName);
210             dictionary newPatch(dictionary::null);
211             newPatch.add("type", boundaryType);
212             newPatch.add("value", ("uniform " + boundaryValue).c_str());
214             oldPatch = newPatch;
215         }
216     }
218     Info<< "    writing updated " << dict.name() << nl << endl;
219     dict.regIOobject::write();
223 void updateCompressibleCase(const fvMesh& mesh)
225     Info<< "Case treated as compressible" << nl << endl;
226     createVolScalarField
227     (
228         mesh,
229         "mut",
230         dimArea/dimTime*dimDensity
231     );
232     replaceBoundaryType
233     (
234         mesh,
235         "mut",
236         compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
237         "0"
238     );
239     replaceBoundaryType
240     (
241         mesh,
242         "epsilon",
243         compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
244             typeName,
245         "0"
246     );
247     replaceBoundaryType
248     (
249         mesh,
250         "omega",
251         compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
252         "0"
253     );
254     replaceBoundaryType
255     (
256         mesh,
257         "k",
258         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
259         "0"
260     );
261     replaceBoundaryType
262     (
263         mesh,
264         "q",
265         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
266         "0"
267     );
268     replaceBoundaryType
269     (
270         mesh,
271         "R",
272         compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
273             typeName,
274         "(0 0 0 0 0 0)"
275     );
279 void updateIncompressibleCase(const fvMesh& mesh)
281     Info<< "Case treated as incompressible" << nl << endl;
282     createVolScalarField(mesh, "nut", dimArea/dimTime);
284     replaceBoundaryType
285     (
286         mesh,
287         "nut",
288         incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
289         "0"
290     );
291     replaceBoundaryType
292     (
293         mesh,
294         "epsilon",
295         incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
296             typeName,
297         "0"
298     );
299     replaceBoundaryType
300     (
301         mesh,
302         "omega",
303         incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
304             typeName,
305         "0"
306     );
307     replaceBoundaryType
308     (
309         mesh,
310         "k",
311         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
312             typeName,
313         "0"
314     );
315     replaceBoundaryType
316     (
317         mesh,
318         "q",
319         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
320             typeName,
321         "0"
322     );
323     replaceBoundaryType
324     (
325         mesh,
326         "R",
327         incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
328             typeName,
329         "(0 0 0 0 0 0)"
330     );
334 int main(int argc, char *argv[])
336     #include "addTimeOptions.H"
337     argList::validOptions.insert("compressible", "");
339     #include "setRootCase.H"
340     #include "createTime.H"
341     #include "createMesh.H"
343     bool compressible = args.optionFound("compressible");
345     Info<< "Updating turbulence fields to operate using new run time "
346         << "selectable" << nl << "wall functions"
347         << nl << endl;
349     if (compressible || caseIsCompressible(mesh))
350     {
351         updateCompressibleCase(mesh);
352     }
353     else
354     {
355         updateIncompressibleCase(mesh);
356     }
358     Info<< "End\n" << endl;
360     return 0;
364 // ************************************************************************* //