Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / preProcessing / applyWallFunctionBoundaryConditions / applyWallFunctionBoundaryConditions.C
blob7710c887c06ce7b5891dbc00cb5d7633366ad6cd
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     applyWallFunctionBounaryConditions
27 Description
28     Updates FOAM 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 "foamTime.H"
38 #include "volFields.H"
39 #include "surfaceFields.H"
41 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
42 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
43 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions/nutWallFunction/nutWallFunctionFvPatchScalarField.H"
44 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
46 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
47 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
48 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions/mutWallFunction/mutWallFunctionFvPatchScalarField.H"
49 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
51 using namespace Foam;
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 bool caseIsCompressible(const fvMesh& mesh)
57     // Attempt flux field
58     IOobject phiHeader
59     (
60         "phi",
61         mesh.time().timeName(),
62         mesh,
63         IOobject::MUST_READ,
64         IOobject::NO_WRITE
65     );
67     if (phiHeader.headerOk())
68     {
69         surfaceScalarField phi(phiHeader, mesh);
70         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
71         {
72             return true;
73         }
74     }
76     // Attempt density field
77     IOobject rhoHeader
78     (
79         "rho",
80         mesh.time().timeName(),
81         mesh,
82         IOobject::MUST_READ,
83         IOobject::NO_WRITE
84     );
86     if (rhoHeader.headerOk())
87     {
88         volScalarField rho(rhoHeader, mesh);
89         if (rho.dimensions() == dimDensity)
90         {
91             return true;
92         }
93     }
95     // Attempt pressure field
96     IOobject pHeader
97     (
98         "p",
99         mesh.time().timeName(),
100         mesh,
101         IOobject::MUST_READ,
102         IOobject::NO_WRITE
103     );
105     if (pHeader.headerOk())
106     {
107         volScalarField p(pHeader, mesh);
108         if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
109         {
110             return true;
111         }
112     }
114     // If none of the above are true, assume that the case is incompressible
115     return false;
119 void createVolScalarField
121     const fvMesh& mesh,
122     const word& fieldName,
123     const dimensionSet& dims
126     IOobject fieldHeader
127     (
128         fieldName,
129         mesh.time().timeName(),
130         mesh,
131         IOobject::MUST_READ,
132         IOobject::NO_WRITE
133     );
135     if (!fieldHeader.headerOk())
136     {
137         Info<< "Creating field " << fieldName << nl << endl;
139         volScalarField field
140         (
141             IOobject
142             (
143                 fieldName,
144                 mesh.time().timeName(),
145                 mesh,
146                 IOobject::NO_READ,
147                 IOobject::NO_WRITE
148             ),
149             mesh,
150             dimensionedScalar("zero", dims, 0.0)
151         );
153         field.write();
154     }
158 void replaceBoundaryType
160     const fvMesh& mesh,
161     const word& fieldName,
162     const word& boundaryType,
163     const string& boundaryValue
166     IOobject header
167     (
168         fieldName,
169         mesh.time().timeName(),
170         mesh,
171         IOobject::MUST_READ,
172         IOobject::NO_WRITE
173     );
175     if (!header.headerOk())
176     {
177         return;
178     }
180     Info<< "Updating boundary types for field " << header.name() << endl;
182     const word oldTypeName = IOdictionary::typeName;
183     const_cast<word&>(IOdictionary::typeName) = word::null;
185     IOdictionary dict(header);
187     const_cast<word&>(IOdictionary::typeName) = oldTypeName;
188     const_cast<word&>(dict.type()) = dict.headerClassName();
190     // Make a backup of the old field
191     word backupName(dict.name() + ".old");
192     Info<< "    copying " << dict.name() << " to "
193         << backupName << endl;
194     IOdictionary dictOld = dict;
195     dictOld.rename(backupName);
196     dictOld.regIOobject::write();
198     // Loop through boundary patches and update
199     const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
200     dictionary& boundaryDict = dict.subDict("boundaryField");
201     forAll(bMesh, patchI)
202     {
203         if (bMesh[patchI].isWall())
204         {
205             word patchName = bMesh[patchI].name();
206             dictionary& oldPatch = boundaryDict.subDict(patchName);
208             dictionary newPatch(dictionary::null);
209             newPatch.add("type", boundaryType);
210             newPatch.add("value", ("uniform " + boundaryValue).c_str());
212             oldPatch = newPatch;
213         }
214     }
216     Info<< "    writing updated " << dict.name() << nl << endl;
217     dict.regIOobject::write();
221 void updateCompressibleCase(const fvMesh& mesh)
223     Info<< "Case treated as compressible" << nl << endl;
224     createVolScalarField
225     (
226         mesh,
227         "mut",
228         dimArea/dimTime*dimDensity
229     );
230     replaceBoundaryType
231     (
232         mesh,
233         "mut",
234         compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
235         "0"
236     );
237     replaceBoundaryType
238     (
239         mesh,
240         "epsilon",
241         compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
242             typeName,
243         "0"
244     );
245     replaceBoundaryType
246     (
247         mesh,
248         "omega",
249         compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
250         "0"
251     );
252     replaceBoundaryType
253     (
254         mesh,
255         "k",
256         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
257         "0"
258     );
259     replaceBoundaryType
260     (
261         mesh,
262         "q",
263         compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
264         "0"
265     );
266     replaceBoundaryType
267     (
268         mesh,
269         "R",
270         compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
271             typeName,
272         "(0 0 0 0 0 0)"
273     );
277 void updateIncompressibleCase(const fvMesh& mesh)
279     Info<< "Case treated as incompressible" << nl << endl;
280     createVolScalarField(mesh, "nut", dimArea/dimTime);
282     replaceBoundaryType
283     (
284         mesh,
285         "nut",
286         incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
287         "0"
288     );
289     replaceBoundaryType
290     (
291         mesh,
292         "epsilon",
293         incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
294             typeName,
295         "0"
296     );
297     replaceBoundaryType
298     (
299         mesh,
300         "omega",
301         incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
302             typeName,
303         "0"
304     );
305     replaceBoundaryType
306     (
307         mesh,
308         "k",
309         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
310             typeName,
311         "0"
312     );
313     replaceBoundaryType
314     (
315         mesh,
316         "q",
317         incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
318             typeName,
319         "0"
320     );
321     replaceBoundaryType
322     (
323         mesh,
324         "R",
325         incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
326             typeName,
327         "(0 0 0 0 0 0)"
328     );
332 int main(int argc, char *argv[])
334     #include "addTimeOptions.H"
335     argList::validOptions.insert("compressible", "");
337     #include "setRootCase.H"
338     #include "createTime.H"
339     #include "createMesh.H"
341     bool compressible = args.optionFound("compressible");
343     Info<< "Updating turbulence fields to operate using new run time "
344         << "selectable" << nl << "wall functions"
345         << nl << endl;
347     if (compressible || caseIsCompressible(mesh))
348     {
349         updateCompressibleCase(mesh);
350     }
351     else
352     {
353         updateIncompressibleCase(mesh);
354     }
356     Info<< "End\n" << endl;
358     return 0;
362 // ************************************************************************* //