1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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/nutkWallFunction/nutkWallFunctionFvPatchScalarField.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/mutkWallFunction/mutkWallFunctionFvPatchScalarField.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 file
193 if (mvBak(dict.objectPath(), "old"))
195 Info<< " Backup original file to "
196 << (dict.objectPath() + ".old") << endl;
199 // Loop through boundary patches and update
200 const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
201 dictionary& boundaryDict = dict.subDict("boundaryField");
202 forAll(bMesh, patchI)
204 if (isA<wallPolyPatch>(bMesh[patchI]))
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());
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;
229 dimArea/dimTime*dimDensity
235 compressible::RASModels::mutkWallFunctionFvPatchScalarField::typeName,
242 compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
250 compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
257 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
264 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
271 compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
278 void updateIncompressibleCase(const fvMesh& mesh)
280 Info<< "Case treated as incompressible" << nl << endl;
281 createVolScalarField(mesh, "nut", dimArea/dimTime);
287 incompressible::RASModels::nutkWallFunctionFvPatchScalarField::typeName,
294 incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
302 incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
310 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
318 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
326 incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
333 int main(int argc, char *argv[])
335 #include "addTimeOptions.H"
336 argList::addBoolOption
339 "force use of compressible wall functions. default is auto-detect."
342 #include "setRootCase.H"
343 #include "createTime.H"
344 #include "createMesh.H"
346 const bool compressible = args.optionFound("compressible");
348 Info<< "Updating turbulence fields to operate using new run time "
349 << "selectable" << nl << "wall functions"
352 if (compressible || caseIsCompressible(mesh))
354 updateCompressibleCase(mesh);
358 updateIncompressibleCase(mesh);
361 Info<< "End\n" << endl;
367 // ************************************************************************* //