FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

backwardsCompatibilityWallFunctions.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "backwardsCompatibilityWallFunctions.H"
00027 
00028 #include <finiteVolume/calculatedFvPatchField.H>
00029 #include <compressibleRASModels/alphatWallFunctionFvPatchScalarField.H>
00030 #include <compressibleRASModels/mutWallFunctionFvPatchScalarField.H>
00031 #include <compressibleRASModels/mutLowReWallFunctionFvPatchScalarField.H>
00032 #include <compressibleRASModels/epsilonWallFunctionFvPatchScalarField.H>
00033 #include <compressibleRASModels/kqRWallFunctionFvPatchField.H>
00034 #include <compressibleRASModels/omegaWallFunctionFvPatchScalarField.H>
00035 
00036 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00037 
00038 namespace Foam
00039 {
00040 namespace compressible
00041 {
00042 
00043 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00044 
00045 tmp<volScalarField> autoCreateAlphat
00046 (
00047     const word& fieldName,
00048     const fvMesh& mesh
00049 )
00050 {
00051     IOobject alphatHeader
00052     (
00053         fieldName,
00054         mesh.time().timeName(),
00055         mesh,
00056         IOobject::MUST_READ,
00057         IOobject::NO_WRITE,
00058         false
00059     );
00060 
00061     if (alphatHeader.headerOk())
00062     {
00063         return tmp<volScalarField>(new volScalarField(alphatHeader, mesh));
00064     }
00065     else
00066     {
00067         Info<< "--> Creating " << fieldName
00068             << " to employ run-time selectable wall functions" << endl;
00069 
00070         const fvBoundaryMesh& bm = mesh.boundary();
00071 
00072         wordList alphatBoundaryTypes(bm.size());
00073 
00074         forAll(bm, patchI)
00075         {
00076             if (isA<wallFvPatch>(bm[patchI]))
00077             {
00078                 alphatBoundaryTypes[patchI] =
00079                     RASModels::alphatWallFunctionFvPatchScalarField::typeName;
00080             }
00081             else
00082             {
00083                 alphatBoundaryTypes[patchI] =
00084                     calculatedFvPatchField<scalar>::typeName;
00085             }
00086         }
00087 
00088         tmp<volScalarField> alphat
00089         (
00090             new volScalarField
00091             (
00092                 IOobject
00093                 (
00094                     fieldName,
00095                     mesh.time().timeName(),
00096                     mesh,
00097                     IOobject::NO_READ,
00098                     IOobject::NO_WRITE,
00099                     false
00100                 ),
00101                 mesh,
00102                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
00103                 alphatBoundaryTypes
00104             )
00105         );
00106 
00107         Info<< "    Writing new " << fieldName << endl;
00108         alphat().write();
00109 
00110         return alphat;
00111     }
00112 }
00113 
00114 
00115 tmp<volScalarField> autoCreateMut
00116 (
00117     const word& fieldName,
00118     const fvMesh& mesh
00119 )
00120 {
00121     IOobject mutHeader
00122     (
00123         fieldName,
00124         mesh.time().timeName(),
00125         mesh,
00126         IOobject::MUST_READ,
00127         IOobject::NO_WRITE,
00128         false
00129     );
00130 
00131     if (mutHeader.headerOk())
00132     {
00133         return tmp<volScalarField>(new volScalarField(mutHeader, mesh));
00134     }
00135     else
00136     {
00137         Info<< "--> Creating " << fieldName
00138             << " to employ run-time selectable wall functions" << endl;
00139 
00140         const fvBoundaryMesh& bm = mesh.boundary();
00141 
00142         wordList mutBoundaryTypes(bm.size());
00143 
00144         forAll(bm, patchI)
00145         {
00146             if (isA<wallFvPatch>(bm[patchI]))
00147             {
00148                 mutBoundaryTypes[patchI] =
00149                     RASModels::mutWallFunctionFvPatchScalarField::typeName;
00150             }
00151             else
00152             {
00153                 mutBoundaryTypes[patchI] =
00154                     calculatedFvPatchField<scalar>::typeName;
00155             }
00156         }
00157 
00158         tmp<volScalarField> mut
00159         (
00160             new volScalarField
00161             (
00162                 IOobject
00163                 (
00164                     fieldName,
00165                     mesh.time().timeName(),
00166                     mesh,
00167                     IOobject::NO_READ,
00168                     IOobject::NO_WRITE,
00169                     false
00170                 ),
00171                 mesh,
00172                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
00173                 mutBoundaryTypes
00174             )
00175         );
00176 
00177         Info<< "    Writing new " << fieldName << endl;
00178         mut().write();
00179 
00180         return mut;
00181     }
00182 }
00183 
00184 
00185 tmp<volScalarField> autoCreateLowReMut
00186 (
00187     const word& fieldName,
00188     const fvMesh& mesh
00189 )
00190 {
00191     IOobject mutHeader
00192     (
00193         fieldName,
00194         mesh.time().timeName(),
00195         mesh,
00196         IOobject::MUST_READ,
00197         IOobject::NO_WRITE,
00198         false
00199     );
00200 
00201     if (mutHeader.headerOk())
00202     {
00203         return tmp<volScalarField>(new volScalarField(mutHeader, mesh));
00204     }
00205     else
00206     {
00207         Info<< "--> Creating " << fieldName
00208             << " to employ run-time selectable wall functions" << endl;
00209 
00210         const fvBoundaryMesh& bm = mesh.boundary();
00211 
00212         wordList mutBoundaryTypes(bm.size());
00213 
00214         forAll(bm, patchI)
00215         {
00216             if (isA<wallFvPatch>(bm[patchI]))
00217             {
00218                 mutBoundaryTypes[patchI] =
00219                     RASModels::mutLowReWallFunctionFvPatchScalarField::typeName;
00220             }
00221             else
00222             {
00223                 mutBoundaryTypes[patchI] =
00224                     calculatedFvPatchField<scalar>::typeName;
00225             }
00226         }
00227 
00228         tmp<volScalarField> mut
00229         (
00230             new volScalarField
00231             (
00232                 IOobject
00233                 (
00234                     fieldName,
00235                     mesh.time().timeName(),
00236                     mesh,
00237                     IOobject::NO_READ,
00238                     IOobject::NO_WRITE,
00239                     false
00240                 ),
00241                 mesh,
00242                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
00243                 mutBoundaryTypes
00244             )
00245         );
00246 
00247         Info<< "    Writing new " << fieldName << endl;
00248         mut().write();
00249 
00250         return mut;
00251     }
00252 }
00253 
00254 
00255 tmp<volScalarField> autoCreateEpsilon
00256 (
00257     const word& fieldName,
00258     const fvMesh& mesh
00259 )
00260 {
00261     return
00262         autoCreateWallFunctionField
00263         <
00264             scalar,
00265             RASModels::epsilonWallFunctionFvPatchScalarField
00266         >
00267         (
00268             fieldName,
00269             mesh
00270         );
00271 }
00272 
00273 
00274 tmp<volScalarField> autoCreateOmega
00275 (
00276     const word& fieldName,
00277     const fvMesh& mesh
00278 )
00279 {
00280     return
00281         autoCreateWallFunctionField
00282         <
00283             scalar,
00284             RASModels::omegaWallFunctionFvPatchScalarField
00285         >
00286         (
00287             fieldName,
00288             mesh
00289         );
00290 }
00291 
00292 
00293 tmp<volScalarField> autoCreateK
00294 (
00295     const word& fieldName,
00296     const fvMesh& mesh
00297 )
00298 {
00299     return
00300         autoCreateWallFunctionField
00301         <
00302             scalar,
00303             RASModels::kqRWallFunctionFvPatchField<scalar>
00304         >
00305         (
00306             fieldName,
00307             mesh
00308         );
00309 }
00310 
00311 
00312 tmp<volScalarField> autoCreateQ
00313 (
00314     const word& fieldName,
00315     const fvMesh& mesh
00316 )
00317 {
00318     return
00319         autoCreateWallFunctionField
00320         <
00321             scalar,
00322             RASModels::kqRWallFunctionFvPatchField<scalar>
00323         >
00324         (
00325             fieldName,
00326             mesh
00327         );
00328 }
00329 
00330 
00331 tmp<volSymmTensorField> autoCreateR
00332 (
00333     const word& fieldName,
00334     const fvMesh& mesh
00335 )
00336 {
00337     return
00338         autoCreateWallFunctionField
00339         <
00340             symmTensor,
00341             RASModels::kqRWallFunctionFvPatchField<symmTensor>
00342         >
00343         (
00344             fieldName,
00345             mesh
00346         );
00347 }
00348 
00349 
00350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00351 
00352 } // End namespace compressible
00353 } // End namespace Foam
00354 
00355 // ************************ vim: set sw=4 sts=4 et: ************************ //
00356 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines