changed isA<wallPolyPatch> & isA<wallFvPatch> almost everywhere
[openfoam-extend-OpenFOAM-1.6-ext.git] / applications / utilities / preProcessing / applyWallFunctionBoundaryConditions / applyWallFunctionBoundaryConditions.C
blobcada9f588d261e217f1a62655c170c388d7c7d20
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
43 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
44 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions/nutWallFunction/nutWallFunctionFvPatchScalarField.H"
45 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
47 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
48 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
49 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions/mutWallFunction/mutWallFunctionFvPatchScalarField.H"
50 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
52 using namespace Foam;
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 bool caseIsCompressible(const fvMesh& mesh)
58     // Attempt flux field
59     IOobject phiHeader
60     (
61         "phi",
62         mesh.time().timeName(),
63         mesh,
64         IOobject::MUST_READ,
65         IOobject::NO_WRITE
66     );
68     if (phiHeader.headerOk())
69     {
70         surfaceScalarField phi(phiHeader, mesh);
71         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
72         {
73             return true;
74         }
75     }
77     // Attempt density field
78     IOobject rhoHeader
79     (
80         "rho",
81         mesh.time().timeName(),
82         mesh,
83         IOobject::MUST_READ,
84         IOobject::NO_WRITE
85     );
87     if (rhoHeader.headerOk())
88     {
89         volScalarField rho(rhoHeader, mesh);
90         if (rho.dimensions() == dimDensity)
91         {
92             return true;
93         }
94     }
96     // Attempt pressure field
97     IOobject pHeader
98     (
99         "p",
100         mesh.time().timeName(),
101         mesh,
102         IOobject::MUST_READ,
103         IOobject::NO_WRITE
104     );
106     if (pHeader.headerOk())
107     {
108         volScalarField p(pHeader, mesh);
109         if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
110         {
111             return true;
112         }
113     }
115     // If none of the above are true, assume that the case is incompressible
116     return false;
120 void createVolScalarField
122     const fvMesh& mesh,
123     const word& fieldName,
124     const dimensionSet& dims
127     IOobject fieldHeader
128     (
129         fieldName,
130         mesh.time().timeName(),
131         mesh,
132         IOobject::MUST_READ,
133         IOobject::NO_WRITE
134     );
136     if (!fieldHeader.headerOk())
137     {
138         Info<< "Creating field " << fieldName << nl << endl;
140         volScalarField field
141         (
142             IOobject
143             (
144                 fieldName,
145                 mesh.time().timeName(),
146                 mesh,
147                 IOobject::NO_READ,
148                 IOobject::NO_WRITE
149             ),
150             mesh,
151             dimensionedScalar("zero", dims, 0.0)
152         );
154         field.write();
155     }
159 void replaceBoundaryType
161     const fvMesh& mesh,
162     const word& fieldName,
163     const word& boundaryType,
164     const string& boundaryValue
167     IOobject header
168     (
169         fieldName,
170         mesh.time().timeName(),
171         mesh,
172         IOobject::MUST_READ,
173         IOobject::NO_WRITE
174     );
176     if (!header.headerOk())
177     {
178         return;
179     }
181     Info<< "Updating boundary types for field " << header.name() << endl;
183     const word oldTypeName = IOdictionary::typeName;
184     const_cast<word&>(IOdictionary::typeName) = word::null;
186     IOdictionary dict(header);
188     const_cast<word&>(IOdictionary::typeName) = oldTypeName;
189     const_cast<word&>(dict.type()) = dict.headerClassName();
191     // Make a backup of the old field
192     word backupName(dict.name() + ".old");
193     Info<< "    copying " << dict.name() << " to "
194         << backupName << endl;
195     IOdictionary dictOld = dict;
196     dictOld.rename(backupName);
197     dictOld.regIOobject::write();
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 (bMesh[patchI].isWall())
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::mutWallFunctionFvPatchScalarField::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::nutWallFunctionFvPatchScalarField::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::validOptions.insert("compressible", "");
338     #include "setRootCase.H"
339     #include "createTime.H"
340     #include "createMesh.H"
342     bool compressible = args.optionFound("compressible");
344     Info<< "Updating turbulence fields to operate using new run time "
345         << "selectable" << nl << "wall functions"
346         << nl << endl;
348     if (compressible || caseIsCompressible(mesh))
349     {
350         updateCompressibleCase(mesh);
351     }
352     else
353     {
354         updateIncompressibleCase(mesh);
355     }
357     Info<< "End\n" << endl;
359     return 0;
363 // ************************************************************************* //