00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068 #include <OpenFOAM/argList.H>
00069 #include <finiteVolume/fvMesh.H>
00070 #include <OpenFOAM/Time.H>
00071 #include <finiteVolume/volFields.H>
00072 #include <finiteVolume/surfaceFields.H>
00073
00074 #include <OpenFOAM/wallPolyPatch.H>
00075
00076 #include <incompressibleRASModels/epsilonWallFunctionFvPatchScalarField.H>
00077 #include <incompressibleRASModels/kqRWallFunctionFvPatchField.H>
00078 #include <incompressibleRASModels/nutWallFunctionFvPatchScalarField.H>
00079 #include <incompressibleRASModels/omegaWallFunctionFvPatchScalarField.H>
00080
00081 #include <compressibleRASModels/epsilonWallFunctionFvPatchScalarField.H>
00082 #include <compressibleRASModels/kqRWallFunctionFvPatchField.H>
00083 #include <compressibleRASModels/mutWallFunctionFvPatchScalarField.H>
00084 #include <compressibleRASModels/omegaWallFunctionFvPatchScalarField.H>
00085
00086 using namespace Foam;
00087
00088
00089
00090 bool caseIsCompressible(const fvMesh& mesh)
00091 {
00092
00093 IOobject phiHeader
00094 (
00095 "phi",
00096 mesh.time().timeName(),
00097 mesh,
00098 IOobject::MUST_READ,
00099 IOobject::NO_WRITE
00100 );
00101
00102 if (phiHeader.headerOk())
00103 {
00104 surfaceScalarField phi(phiHeader, mesh);
00105 if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
00106 {
00107 return true;
00108 }
00109 }
00110
00111
00112 IOobject rhoHeader
00113 (
00114 "rho",
00115 mesh.time().timeName(),
00116 mesh,
00117 IOobject::MUST_READ,
00118 IOobject::NO_WRITE
00119 );
00120
00121 if (rhoHeader.headerOk())
00122 {
00123 volScalarField rho(rhoHeader, mesh);
00124 if (rho.dimensions() == dimDensity)
00125 {
00126 return true;
00127 }
00128 }
00129
00130
00131 IOobject pHeader
00132 (
00133 "p",
00134 mesh.time().timeName(),
00135 mesh,
00136 IOobject::MUST_READ,
00137 IOobject::NO_WRITE
00138 );
00139
00140 if (pHeader.headerOk())
00141 {
00142 volScalarField p(pHeader, mesh);
00143 if (p.dimensions() == dimMass/sqr(dimTime)/dimLength)
00144 {
00145 return true;
00146 }
00147 }
00148
00149
00150 return false;
00151 }
00152
00153
00154 void createVolScalarField
00155 (
00156 const fvMesh& mesh,
00157 const word& fieldName,
00158 const dimensionSet& dims
00159 )
00160 {
00161 IOobject fieldHeader
00162 (
00163 fieldName,
00164 mesh.time().timeName(),
00165 mesh,
00166 IOobject::MUST_READ,
00167 IOobject::NO_WRITE
00168 );
00169
00170 if (!fieldHeader.headerOk())
00171 {
00172 Info<< "Creating field " << fieldName << nl << endl;
00173
00174 volScalarField field
00175 (
00176 IOobject
00177 (
00178 fieldName,
00179 mesh.time().timeName(),
00180 mesh,
00181 IOobject::NO_READ,
00182 IOobject::NO_WRITE
00183 ),
00184 mesh,
00185 dimensionedScalar("zero", dims, 0.0)
00186 );
00187
00188 field.write();
00189 }
00190 }
00191
00192
00193 void replaceBoundaryType
00194 (
00195 const fvMesh& mesh,
00196 const word& fieldName,
00197 const word& boundaryType,
00198 const string& boundaryValue
00199 )
00200 {
00201 IOobject header
00202 (
00203 fieldName,
00204 mesh.time().timeName(),
00205 mesh,
00206 IOobject::MUST_READ,
00207 IOobject::NO_WRITE
00208 );
00209
00210 if (!header.headerOk())
00211 {
00212 return;
00213 }
00214
00215 Info<< "Updating boundary types for field " << header.name() << endl;
00216
00217 const word oldTypeName = IOdictionary::typeName;
00218 const_cast<word&>(IOdictionary::typeName) = word::null;
00219
00220 IOdictionary dict(header);
00221
00222 const_cast<word&>(IOdictionary::typeName) = oldTypeName;
00223 const_cast<word&>(dict.type()) = dict.headerClassName();
00224
00225
00226 word backupName(dict.name() + ".old");
00227 Info<< " copying " << dict.name() << " to "
00228 << backupName << endl;
00229 IOdictionary dictOld = dict;
00230 dictOld.rename(backupName);
00231 dictOld.regIOobject::write();
00232
00233
00234 const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
00235 dictionary& boundaryDict = dict.subDict("boundaryField");
00236 forAll(bMesh, patchI)
00237 {
00238 if (isA<wallPolyPatch>(bMesh[patchI]))
00239 {
00240 word patchName = bMesh[patchI].name();
00241 dictionary& oldPatch = boundaryDict.subDict(patchName);
00242
00243 dictionary newPatch(dictionary::null);
00244 newPatch.add("type", boundaryType);
00245 newPatch.add("value", ("uniform " + boundaryValue).c_str());
00246
00247 oldPatch = newPatch;
00248 }
00249 }
00250
00251 Info<< " writing updated " << dict.name() << nl << endl;
00252 dict.regIOobject::write();
00253 }
00254
00255
00256 void updateCompressibleCase(const fvMesh& mesh)
00257 {
00258 Info<< "Case treated as compressible" << nl << endl;
00259 createVolScalarField
00260 (
00261 mesh,
00262 "mut",
00263 dimArea/dimTime*dimDensity
00264 );
00265 replaceBoundaryType
00266 (
00267 mesh,
00268 "mut",
00269 compressible::RASModels::mutWallFunctionFvPatchScalarField::typeName,
00270 "0"
00271 );
00272 replaceBoundaryType
00273 (
00274 mesh,
00275 "epsilon",
00276 compressible::RASModels::epsilonWallFunctionFvPatchScalarField::
00277 typeName,
00278 "0"
00279 );
00280 replaceBoundaryType
00281 (
00282 mesh,
00283 "omega",
00284 compressible::RASModels::omegaWallFunctionFvPatchScalarField::typeName,
00285 "0"
00286 );
00287 replaceBoundaryType
00288 (
00289 mesh,
00290 "k",
00291 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
00292 "0"
00293 );
00294 replaceBoundaryType
00295 (
00296 mesh,
00297 "q",
00298 compressible::RASModels::kqRWallFunctionFvPatchField<scalar>::typeName,
00299 "0"
00300 );
00301 replaceBoundaryType
00302 (
00303 mesh,
00304 "R",
00305 compressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
00306 typeName,
00307 "(0 0 0 0 0 0)"
00308 );
00309 }
00310
00311
00312 void updateIncompressibleCase(const fvMesh& mesh)
00313 {
00314 Info<< "Case treated as incompressible" << nl << endl;
00315 createVolScalarField(mesh, "nut", dimArea/dimTime);
00316
00317 replaceBoundaryType
00318 (
00319 mesh,
00320 "nut",
00321 incompressible::RASModels::nutWallFunctionFvPatchScalarField::typeName,
00322 "0"
00323 );
00324 replaceBoundaryType
00325 (
00326 mesh,
00327 "epsilon",
00328 incompressible::RASModels::epsilonWallFunctionFvPatchScalarField::
00329 typeName,
00330 "0"
00331 );
00332 replaceBoundaryType
00333 (
00334 mesh,
00335 "omega",
00336 incompressible::RASModels::omegaWallFunctionFvPatchScalarField::
00337 typeName,
00338 "0"
00339 );
00340 replaceBoundaryType
00341 (
00342 mesh,
00343 "k",
00344 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
00345 typeName,
00346 "0"
00347 );
00348 replaceBoundaryType
00349 (
00350 mesh,
00351 "q",
00352 incompressible::RASModels::kqRWallFunctionFvPatchField<scalar>::
00353 typeName,
00354 "0"
00355 );
00356 replaceBoundaryType
00357 (
00358 mesh,
00359 "R",
00360 incompressible::RASModels::kqRWallFunctionFvPatchField<symmTensor>::
00361 typeName,
00362 "(0 0 0 0 0 0)"
00363 );
00364 }
00365
00366
00367 int main(int argc, char *argv[])
00368 {
00369 #include <OpenFOAM/addTimeOptions.H>
00370 argList::validOptions.insert("compressible", "");
00371
00372 #include <OpenFOAM/setRootCase.H>
00373 #include <OpenFOAM/createTime.H>
00374 #include <OpenFOAM/createMesh.H>
00375
00376 bool compressible = args.optionFound("compressible");
00377
00378 Info<< "Updating turbulence fields to operate using new run time "
00379 << "selectable" << nl << "wall functions"
00380 << nl << endl;
00381
00382 if (compressible || caseIsCompressible(mesh))
00383 {
00384 updateCompressibleCase(mesh);
00385 }
00386 else
00387 {
00388 updateIncompressibleCase(mesh);
00389 }
00390
00391 Info<< "End\n" << endl;
00392
00393 return 0;
00394 }
00395
00396
00397