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

yPlusRAS.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) 1991-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 Application
00025     yPlusRAS
00026 
00027 Description
00028     Calculates and reports yPlus for all wall patches, for the specified times
00029     when using RAS turbulence models.
00030 
00031     Default behaviour assumes operating in incompressible mode. To apply to
00032     compressible RAS cases, use the -compressible option.
00033 
00034 Usage
00035 
00036     - yPlusRAS [OPTIONS]
00037 
00038     @param -compressible
00039     Operate in compressible mode.
00040 
00041     @param -noZero \n
00042     Ignore timestep 0.
00043 
00044     @param -constant \n
00045     Include the constant directory.
00046 
00047     @param -time <time>\n
00048     Apply only to specific time.
00049 
00050     @param -latestTime \n
00051     Only apply to latest time step.
00052 
00053     @param -case <dir>\n
00054     Case directory.
00055 
00056     @param -parallel \n
00057     Run in parallel.
00058 
00059     @param -help \n
00060     Display help message.
00061 
00062     @param -doc \n
00063     Display Doxygen API documentation page for this application.
00064 
00065     @param -srcDoc \n
00066     Display Doxygen source documentation page for this application.
00067 
00068 \*---------------------------------------------------------------------------*/
00069 
00070 #include <finiteVolume/fvCFD.H>
00071 
00072 #include <incompressibleTransportModels/singlePhaseTransportModel.H>
00073 #include <incompressibleRASModels/RASModel.H>
00074 #include <incompressibleRASModels/nutWallFunctionFvPatchScalarField.H>
00075 
00076 #include <basicThermophysicalModels/basicPsiThermo.H>
00077 #include <compressibleRASModels/RASModel.H>
00078 #include <compressibleRASModels/mutWallFunctionFvPatchScalarField.H>
00079 
00080 #include <finiteVolume/wallDist.H>
00081 
00082 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00083 
00084 void calcIncompressibleYPlus
00085 (
00086     const fvMesh& mesh,
00087     const Time& runTime,
00088     const volVectorField& U,
00089     volScalarField& yPlus
00090 )
00091 {
00092     typedef incompressible::RASModels::nutWallFunctionFvPatchScalarField
00093         wallFunctionPatchField;
00094 
00095     #include <finiteVolume/createPhi.H>
00096 
00097     singlePhaseTransportModel laminarTransport(U, phi);
00098 
00099     autoPtr<incompressible::RASModel> RASModel
00100     (
00101         incompressible::RASModel::New(U, phi, laminarTransport)
00102     );
00103 
00104     const volScalarField::GeometricBoundaryField nutPatches =
00105         RASModel->nut()().boundaryField();
00106 
00107     bool foundNutPatch = false;
00108     forAll(nutPatches, patchi)
00109     {
00110         if (isA<wallFunctionPatchField>(nutPatches[patchi]))
00111         {
00112             foundNutPatch = true;
00113 
00114             const wallFunctionPatchField& nutPw =
00115                 dynamic_cast<const wallFunctionPatchField&>
00116                     (nutPatches[patchi]);
00117 
00118             yPlus.boundaryField()[patchi] = nutPw.yPlus();
00119             const scalarField& Yp = yPlus.boundaryField()[patchi];
00120 
00121             Info<< "Patch " << patchi
00122                 << " named " << nutPw.patch().name()
00123                 << " y+ : min: " << min(Yp) << " max: " << max(Yp)
00124                 << " average: " << average(Yp) << nl << endl;
00125         }
00126     }
00127 
00128     if (!foundNutPatch)
00129     {
00130         Info<< "    no " << wallFunctionPatchField::typeName << " patches"
00131             << endl;
00132     }
00133 }
00134 
00135 
00136 void calcCompressibleYPlus
00137 (
00138     const fvMesh& mesh,
00139     const Time& runTime,
00140     const volVectorField& U,
00141     volScalarField& yPlus
00142 )
00143 {
00144     typedef compressible::RASModels::mutWallFunctionFvPatchScalarField
00145         wallFunctionPatchField;
00146 
00147     IOobject rhoHeader
00148     (
00149         "rho",
00150         runTime.timeName(),
00151         mesh,
00152         IOobject::MUST_READ,
00153         IOobject::NO_WRITE
00154     );
00155 
00156     if (!rhoHeader.headerOk())
00157     {
00158         Info<< "    no rho field" << endl;
00159         return;
00160     }
00161 
00162     Info << "Reading field rho\n" << endl;
00163     volScalarField rho(rhoHeader, mesh);
00164 
00165     #include <finiteVolume/compressibleCreatePhi.H>
00166 
00167     autoPtr<basicPsiThermo> pThermo
00168     (
00169         basicPsiThermo::New(mesh)
00170     );
00171     basicPsiThermo& thermo = pThermo();
00172 
00173     autoPtr<compressible::RASModel> RASModel
00174     (
00175         compressible::RASModel::New
00176         (
00177             rho,
00178             U,
00179             phi,
00180             thermo
00181         )
00182     );
00183 
00184     const volScalarField::GeometricBoundaryField mutPatches =
00185         RASModel->mut()().boundaryField();
00186 
00187     bool foundMutPatch = false;
00188     forAll(mutPatches, patchi)
00189     {
00190         if (isA<wallFunctionPatchField>(mutPatches[patchi]))
00191         {
00192             foundMutPatch = true;
00193 
00194             const wallFunctionPatchField& mutPw =
00195                 dynamic_cast<const wallFunctionPatchField&>
00196                     (mutPatches[patchi]);
00197 
00198             yPlus.boundaryField()[patchi] = mutPw.yPlus();
00199             const scalarField& Yp = yPlus.boundaryField()[patchi];
00200 
00201             Info<< "Patch " << patchi
00202                 << " named " << mutPw.patch().name()
00203                 << " y+ : min: " << min(Yp) << " max: " << max(Yp)
00204                 << " average: " << average(Yp) << nl << endl;
00205         }
00206     }
00207 
00208     if (!foundMutPatch)
00209     {
00210         Info<< "    no " << wallFunctionPatchField::typeName << " patches"
00211             << endl;
00212     }
00213 }
00214 
00215 
00216 int main(int argc, char *argv[])
00217 {
00218     timeSelector::addOptions();
00219 
00220     #include <OpenFOAM/addRegionOption.H>
00221 
00222     argList::validOptions.insert("compressible","");
00223 
00224     #include <OpenFOAM/setRootCase.H>
00225     #include <OpenFOAM/createTime.H>
00226     instantList timeDirs = timeSelector::select0(runTime, args);
00227     #include <OpenFOAM/createNamedMesh.H>
00228 
00229     bool compressible = args.optionFound("compressible");
00230 
00231     forAll(timeDirs, timeI)
00232     {
00233         runTime.setTime(timeDirs[timeI], timeI);
00234         Info<< "Time = " << runTime.timeName() << endl;
00235         fvMesh::readUpdateState state = mesh.readUpdate();
00236 
00237         // Wall distance
00238         if (timeI == 0 || state != fvMesh::UNCHANGED)
00239         {
00240             Info<< "Calculating wall distance\n" << endl;
00241             wallDist y(mesh, true);
00242             Info<< "Writing wall distance to field "
00243                 << y.name() << nl << endl;
00244             y.write();
00245         }
00246 
00247         volScalarField yPlus
00248         (
00249             IOobject
00250             (
00251                 "yPlus",
00252                 runTime.timeName(),
00253                 mesh,
00254                 IOobject::NO_READ,
00255                 IOobject::NO_WRITE
00256             ),
00257             mesh,
00258             dimensionedScalar("yPlus", dimless, 0.0)
00259         );
00260 
00261         IOobject UHeader
00262         (
00263             "U",
00264             runTime.timeName(),
00265             mesh,
00266             IOobject::MUST_READ,
00267             IOobject::NO_WRITE
00268         );
00269 
00270         if (UHeader.headerOk())
00271         {
00272             Info << "Reading field U\n" << endl;
00273             volVectorField U(UHeader, mesh);
00274 
00275             if (compressible)
00276             {
00277                 calcCompressibleYPlus(mesh, runTime, U, yPlus);
00278             }
00279             else
00280             {
00281                 calcIncompressibleYPlus(mesh, runTime, U, yPlus);
00282             }
00283         }
00284         else
00285         {
00286             Info<< "    no U field" << endl;
00287         }
00288 
00289         Info<< "Writing yPlus to field " << yPlus.name() << nl << endl;
00290 
00291         yPlus.write();
00292     }
00293 
00294     Info<< "End\n" << endl;
00295 
00296     return 0;
00297 }
00298 
00299 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines