1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
25 applyWallFunctionBounaryConditions
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 \*---------------------------------------------------------------------------*/
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"
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 bool caseIsCompressible(const fvMesh& mesh)
61 mesh.time().timeName(),
67 if (phiHeader.headerOk())
69 surfaceScalarField phi(phiHeader, mesh);
70 if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
76 // Attempt density field
80 mesh.time().timeName(),
86 if (rhoHeader.headerOk())
88 volScalarField rho(rhoHeader, mesh);
89 if (rho.dimensions() == dimDensity)
95 // Attempt pressure field
99 mesh.time().timeName(),
105 if (pHeader.headerOk())
107 volScalarField p(pHeader, mesh);
108 if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
114 // If none of the above are true, assume that the case is incompressible
119 void createVolScalarField
122 const word& fieldName,
123 const dimensionSet& dims
129 mesh.time().timeName(),
135 if (!fieldHeader.headerOk())
137 Info<< "Creating field " << fieldName << nl << endl;
144 mesh.time().timeName(),
150 dimensionedScalar("zero", dims, 0.0)
158 void replaceBoundaryType
161 const word& fieldName,
162 const word& boundaryType,
163 const string& boundaryValue
169 mesh.time().timeName(),
175 if (!header.headerOk())
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)
203 if (bMesh[patchI].isWall())
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());
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;
228 dimArea/dimTime*dimDensity
234 compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
241 compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
249 compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
256 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
263 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
270 compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
277 void updateIncompressibleCase(const fvMesh& mesh)
279 Info<< "Case treated as incompressible" << nl << endl;
280 createVolScalarField(mesh, "nut", dimArea/dimTime);
286 incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
293 incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
301 incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
309 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
317 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
325 incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
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"
347 if (compressible || caseIsCompressible(mesh))
349 updateCompressibleCase(mesh);
353 updateIncompressibleCase(mesh);
356 Info<< "End\n" << endl;
362 // ************************************************************************* //