1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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/>.
25 applyWallFunctionBounaryConditions
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 \*---------------------------------------------------------------------------*/
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"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 bool caseIsCompressible(const fvMesh& mesh)
63 mesh.time().timeName(),
69 if (phiHeader.headerOk())
71 surfaceScalarField phi(phiHeader, mesh);
72 if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
78 // Attempt density field
82 mesh.time().timeName(),
88 if (rhoHeader.headerOk())
90 volScalarField rho(rhoHeader, mesh);
91 if (rho.dimensions() == dimDensity)
97 // Attempt pressure field
101 mesh.time().timeName(),
107 if (pHeader.headerOk())
109 volScalarField p(pHeader, mesh);
110 if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
116 // If none of the above are true, assume that the case is incompressible
121 void createVolScalarField
124 const word& fieldName,
125 const dimensionSet& dims
131 mesh.time().timeName(),
137 if (!fieldHeader.headerOk())
139 Info<< "Creating field " << fieldName << nl << endl;
146 mesh.time().timeName(),
152 dimensionedScalar("zero", dims, 0.0)
160 void replaceBoundaryType
163 const word& fieldName,
164 const word& boundaryType,
165 const string& boundaryValue
171 mesh.time().timeName(),
177 if (!header.headerOk())
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)
205 if (isA<wallPolyPatch>(bMesh[patchI]))
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());
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;
230 dimArea/dimTime*dimDensity
236 compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
243 compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
251 compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
258 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
265 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
272 compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
279 void updateIncompressibleCase(const fvMesh& mesh)
281 Info<< "Case treated as incompressible" << nl << endl;
282 createVolScalarField(mesh, "nut", dimArea/dimTime);
288 incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
295 incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
303 incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
311 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
319 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
327 incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
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"
349 if (compressible || caseIsCompressible(mesh))
351 updateCompressibleCase(mesh);
355 updateIncompressibleCase(mesh);
358 Info<< "End\n" << endl;
364 // ************************************************************************* //