initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / preProcessing / applyWallFunctionBoundaryConditions / applyWallFunctionBoundaryConditions.C
blob166386c8e3333bbf4c61e8928d982316e4df5290
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     applyWallFunctionBounaryConditions
28 Description
29     Updates OpenFOAM RAS cases to use the new (v1.6) wall function framework
31     Attempts to determine whether case is compressible or incompressible, or
32     can be supplied with -compressible command line argument
34 \*---------------------------------------------------------------------------*/
36 #include "argList.H"
37 #include "fvMesh.H"
38 #include "Time.H"
39 #include "volFields.H"
40 #include "surfaceFields.H"
42 #include "wallPolyPatch.H"
44 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
45 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
46 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions/nutWallFunction/nutWallFunctionFvPatchScalarField.H"
47 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
49 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
50 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
51 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions/mutWallFunction/mutWallFunctionFvPatchScalarField.H"
52 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
54 using namespace Foam;
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 bool caseIsCompressible(const fvMesh& mesh)
60     // Attempt flux field
61     IOobject phiHeader
62     (
63         "phi",
64         mesh.time().timeName(),
65         mesh,
66         IOobject::MUST_READ,
67         IOobject::NO_WRITE
68     );
70     if (phiHeader.headerOk())
71     {
72         surfaceScalarField phi(phiHeader, mesh);
73         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
74         {
75             return true;
76         }
77     }
79     // Attempt density field
80     IOobject rhoHeader
81     (
82         "rho",
83         mesh.time().timeName(),
84         mesh,
85         IOobject::MUST_READ,
86         IOobject::NO_WRITE
87     );
89     if (rhoHeader.headerOk())
90     {
91         volScalarField rho(rhoHeader, mesh);
92         if (rho.dimensions() == dimDensity)
93         {
94             return true;
95         }
96     }
98     // Attempt pressure field
99     IOobject pHeader
100     (
101         "p",
102         mesh.time().timeName(),
103         mesh,
104         IOobject::MUST_READ,
105         IOobject::NO_WRITE
106     );
108     if (pHeader.headerOk())
109     {
110         volScalarField p(pHeader, mesh);
111         if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
112         {
113             return true;
114         }
115     }
117     // If none of the above are true, assume that the case is incompressible
118     return false;
122 void createVolScalarField
124     const fvMesh& mesh,
125     const word& fieldName,
126     const dimensionSet& dims
129     IOobject fieldHeader
130     (
131         fieldName,
132         mesh.time().timeName(),
133         mesh,
134         IOobject::MUST_READ,
135         IOobject::NO_WRITE
136     );
138     if (!fieldHeader.headerOk())
139     {
140         Info<< "Creating field " << fieldName << nl << endl;
142         volScalarField field
143         (
144             IOobject
145             (
146                 fieldName,
147                 mesh.time().timeName(),
148                 mesh,
149                 IOobject::NO_READ,
150                 IOobject::NO_WRITE
151             ),
152             mesh,
153             dimensionedScalar("zero", dims, 0.0)
154         );
156         field.write();
157     }
161 void replaceBoundaryType
163     const fvMesh& mesh,
164     const word& fieldName,
165     const word& boundaryType,
166     const string& boundaryValue
169     IOobject header
170     (
171         fieldName,
172         mesh.time().timeName(),
173         mesh,
174         IOobject::MUST_READ,
175         IOobject::NO_WRITE
176     );
178     if (!header.headerOk())
179     {
180         return;
181     }
183     Info<< "Updating boundary types for field " << header.name() << endl;
185     const word oldTypeName = IOdictionary::typeName;
186     const_cast<word&>(IOdictionary::typeName) = word::null;
188     IOdictionary dict(header);
190     const_cast<word&>(IOdictionary::typeName) = oldTypeName;
191     const_cast<word&>(dict.type()) = dict.headerClassName();
193     // Make a backup of the old field
194     word backupName(dict.name() + ".old");
195     Info<< "    copying " << dict.name() << " to "
196         << backupName << endl;
197     IOdictionary dictOld = dict;
198     dictOld.rename(backupName);
199     dictOld.regIOobject::write();
201     // Loop through boundary patches and update
202     const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
203     dictionary& boundaryDict = dict.subDict("boundaryField");
204     forAll(bMesh, patchI)
205     {
206         if (isType<wallPolyPatch>(bMesh[patchI]))
207         {
208             word patchName = bMesh[patchI].name();
209             dictionary& oldPatch = boundaryDict.subDict(patchName);
211             dictionary newPatch(dictionary::null);
212             newPatch.add("type", boundaryType);
213             newPatch.add("value", ("uniform " + boundaryValue).c_str());
215             oldPatch = newPatch;
216         }
217     }
219     Info<< "    writing updated " << dict.name() << nl << endl;
220     dict.regIOobject::write();
224 void updateCompressibleCase(const fvMesh& mesh)
226     Info<< "Case treated as compressible" << nl << endl;
227     createVolScalarField
228     (
229         mesh,
230         "mut",
231         dimArea/dimTime*dimDensity
232     );
233     replaceBoundaryType
234     (
235         mesh,
236         "mut",
237         compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
238         "0"
239     );
240     replaceBoundaryType
241     (
242         mesh,
243         "epsilon",
244         compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
245             typeName,
246         "0"
247     );
248     replaceBoundaryType
249     (
250         mesh,
251         "omega",
252         compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
253         "0"
254     );
255     replaceBoundaryType
256     (
257         mesh,
258         "k",
259         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
260         "0"
261     );
262     replaceBoundaryType
263     (
264         mesh,
265         "q",
266         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
267         "0"
268     );
269     replaceBoundaryType
270     (
271         mesh,
272         "R",
273         compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
274             typeName,
275         "(0 0 0 0 0 0)"
276     );
280 void updateIncompressibleCase(const fvMesh& mesh)
282     Info<< "Case treated as incompressible" << nl << endl;
283     createVolScalarField(mesh, "nut", dimArea/dimTime);
285     replaceBoundaryType
286     (
287         mesh,
288         "nut",
289         incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
290         "0"
291     );
292     replaceBoundaryType
293     (
294         mesh,
295         "epsilon",
296         incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
297             typeName,
298         "0"
299     );
300     replaceBoundaryType
301     (
302         mesh,
303         "omega",
304         incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
305             typeName,
306         "0"
307     );
308     replaceBoundaryType
309     (
310         mesh,
311         "k",
312         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
313             typeName,
314         "0"
315     );
316     replaceBoundaryType
317     (
318         mesh,
319         "q",
320         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
321             typeName,
322         "0"
323     );
324     replaceBoundaryType
325     (
326         mesh,
327         "R",
328         incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
329             typeName,
330         "(0 0 0 0 0 0)"
331     );
335 int main(int argc, char *argv[])
337     #include "addTimeOptions.H"
338     argList::validOptions.insert("compressible", "");
340     #include "setRootCase.H"
341     #include "createTime.H"
342     #include "createMesh.H"
344     bool compressible = args.optionFound("compressible");
346     Info<< "Updating turbulence fields to operate using new run time "
347         << "selectable" << nl << "wall functions"
348         << nl << endl;
350     if (compressible || caseIsCompressible(mesh))
351     {
352         updateCompressibleCase(mesh);
353     }
354     else
355     {
356         updateIncompressibleCase(mesh);
357     }
359     Info<< "End\n" << endl;
361     return 0;
365 // ************************************************************************* //