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 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
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
26 applyWallFunctionBounaryConditions
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 \*---------------------------------------------------------------------------*/
39 #include "volFields.H"
40 #include "surfaceFields.H"
42 #include "wallPolyPatch.H"
44 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
45 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
46 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions/nutWallFunction/nutWallFunctionFvPatchScalarField.H"
47 #include "incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
49 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.H"
50 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/kqRWallFunctions/kqRWallFunction/kqRWallFunctionFvPatchField.H"
51 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions/mutWallFunction/mutWallFunctionFvPatchScalarField.H"
52 #include "compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.H"
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 bool caseIsCompressible(const fvMesh& mesh)
64 mesh.time().timeName(),
70 if (phiHeader.headerOk())
72 surfaceScalarField phi(phiHeader, mesh);
73 if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
79 // Attempt density field
83 mesh.time().timeName(),
89 if (rhoHeader.headerOk())
91 volScalarField rho(rhoHeader, mesh);
92 if (rho.dimensions() == dimDensity)
98 // Attempt pressure field
102 mesh.time().timeName(),
108 if (pHeader.headerOk())
110 volScalarField p(pHeader, mesh);
111 if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
117 // If none of the above are true, assume that the case is incompressible
122 void createVolScalarField
125 const word& fieldName,
126 const dimensionSet& dims
132 mesh.time().timeName(),
138 if (!fieldHeader.headerOk())
140 Info<< "Creating field " << fieldName << nl << endl;
147 mesh.time().timeName(),
153 dimensionedScalar("zero", dims, 0.0)
161 void replaceBoundaryType
164 const word& fieldName,
165 const word& boundaryType,
166 const string& boundaryValue
172 mesh.time().timeName(),
178 if (!header.headerOk())
183 Info<< "Updating boundary types for field " << header.name() << endl;
185 const word oldTypeName = IOdictionary::typeName;
186 const_cast<word&>(IOdictionary::typeName) = word::null;
188 IOdictionary dict(header);
190 const_cast<word&>(IOdictionary::typeName) = oldTypeName;
191 const_cast<word&>(dict.type()) = dict.headerClassName();
193 // Make a backup of the old field
194 word backupName(dict.name() + ".old");
195 Info<< " copying " << dict.name() << " to "
196 << backupName << endl;
197 IOdictionary dictOld = dict;
198 dictOld.rename(backupName);
199 dictOld.regIOobject::write();
201 // Loop through boundary patches and update
202 const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
203 dictionary& boundaryDict = dict.subDict("boundaryField");
204 forAll(bMesh, patchI)
206 if (isType<wallPolyPatch>(bMesh[patchI]))
208 word patchName = bMesh[patchI].name();
209 dictionary& oldPatch = boundaryDict.subDict(patchName);
211 dictionary newPatch(dictionary::null);
212 newPatch.add("type", boundaryType);
213 newPatch.add("value", ("uniform " + boundaryValue).c_str());
219 Info<< " writing updated " << dict.name() << nl << endl;
220 dict.regIOobject::write();
224 void updateCompressibleCase(const fvMesh& mesh)
226 Info<< "Case treated as compressible" << nl << endl;
231 dimArea/dimTime*dimDensity
237 compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
244 compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
252 compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
259 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
266 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
273 compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
280 void updateIncompressibleCase(const fvMesh& mesh)
282 Info<< "Case treated as incompressible" << nl << endl;
283 createVolScalarField(mesh, "nut", dimArea/dimTime);
289 incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
296 incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
304 incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
312 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
320 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
328 incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
335 int main(int argc, char *argv[])
337 #include "addTimeOptions.H"
338 argList::validOptions.insert("compressible", "");
340 #include "setRootCase.H"
341 #include "createTime.H"
342 #include "createMesh.H"
344 bool compressible = args.optionFound("compressible");
346 Info<< "Updating turbulence fields to operate using new run time "
347 << "selectable" << nl << "wall functions"
350 if (compressible || caseIsCompressible(mesh))
352 updateCompressibleCase(mesh);
356 updateIncompressibleCase(mesh);
359 Info<< "End\n" << endl;
365 // ************************************************************************* //