ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / wall / yPlusLES / yPlusLES.C
blob0bc6773518e380a46abf872e27e8633af0e7efcc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     yPlusLES
27 Description
28     Calculates and reports yPlus for all wall patches, for the specified times.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
33 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
34 #include "LESModel.H"
35 #include "nearWallDist.H"
36 #include "wallDist.H"
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 int main(int argc, char *argv[])
42     timeSelector::addOptions();
43     #include "setRootCase.H"
44 #   include "createTime.H"
45     instantList timeDirs = timeSelector::select0(runTime, args);
46 #   include "createMesh.H"
48     forAll(timeDirs, timeI)
49     {
50         runTime.setTime(timeDirs[timeI], timeI);
51         Info<< "Time = " << runTime.timeName() << endl;
52         fvMesh::readUpdateState state = mesh.readUpdate();
54         // Wall distance
55         if (timeI == 0 || state != fvMesh::UNCHANGED)
56         {
57             Info<< "Calculating wall distance\n" << endl;
58             wallDist y(mesh, true);
59             Info<< "Writing wall distance to field "
60                 << y.name() << nl << endl;
61             y.write();
62         }
65         volScalarField yPlus
66         (
67             IOobject
68             (
69                 "yPlus",
70                 runTime.timeName(),
71                 mesh,
72                 IOobject::NO_READ,
73                 IOobject::NO_WRITE
74             ),
75             mesh,
76             dimensionedScalar("yPlus", dimless, 0.0)
77         );
79         Info<< "Reading field U\n" << endl;
80         volVectorField U
81         (
82             IOobject
83             (
84                 "U",
85                 runTime.timeName(),
86                 mesh,
87                 IOobject::MUST_READ,
88                 IOobject::NO_WRITE
89             ),
90             mesh
91         );
93 #       include "createPhi.H"
95         singlePhaseTransportModel laminarTransport(U, phi);
97         autoPtr<incompressible::LESModel> sgsModel
98         (
99             incompressible::LESModel::New(U, phi, laminarTransport)
100         );
102         volScalarField::GeometricBoundaryField d = nearWallDist(mesh).y();
103         volScalarField nuEff(sgsModel->nuEff());
105         const fvPatchList& patches = mesh.boundary();
107         const volScalarField nuLam(sgsModel->nu());
109         forAll(patches, patchi)
110         {
111             const fvPatch& currPatch = patches[patchi];
113             if (isA<wallFvPatch>(currPatch))
114             {
115                 yPlus.boundaryField()[patchi] =
116                     d[patchi]
117                    *sqrt
118                     (
119                         nuEff.boundaryField()[patchi]
120                        *mag(U.boundaryField()[patchi].snGrad())
121                     )
122                    /nuLam.boundaryField()[patchi];
123                 const scalarField& Yp = yPlus.boundaryField()[patchi];
125                 Info<< "Patch " << patchi
126                     << " named " << currPatch.name()
127                     << " y+ : min: " << min(Yp) << " max: " << max(Yp)
128                     << " average: " << average(Yp) << nl << endl;
129             }
130         }
132         Info<< "Writing yPlus to field "
133             << yPlus.name() << nl << endl;
135         yPlus.write();
136     }
138     Info<< "End\n" << endl;
140     return 0;
144 // ************************************************************************* //