changed isA<wallPolyPatch> & isA<wallFvPatch> almost everywhere
[openfoam-extend-OpenFOAM-1.6-ext.git] / applications / utilities / postProcessing / wall / yPlusLES / yPlusLES.C
blob2c2f1412ce27e7cae3d27e86ff488b7fd59a53cd
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 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
19     for more details.
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
25 Application
26     yPlusLES
28 Description
29     Calculates and reports yPlus for all wall patches, for the specified times.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
35 #include "LESModel.H"
36 #include "nearWallDist.H"
37 #include "wallDist.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     timeSelector::addOptions();
44     #include "setRootCase.H"
45 #   include "createTime.H"
46     instantList timeDirs = timeSelector::select0(runTime, args);
47 #   include "createMesh.H"
49     forAll(timeDirs, timeI)
50     {
51         runTime.setTime(timeDirs[timeI], timeI);
52         Info<< "Time = " << runTime.timeName() << endl;
53         fvMesh::readUpdateState state = mesh.readUpdate();
55         // Wall distance
56         if (timeI == 0 || state != fvMesh::UNCHANGED)
57         {
58             Info<< "Calculating wall distance\n" << endl;
59             wallDist y(mesh, true);
60             Info<< "Writing wall distance to field "
61                 << y.name() << nl << endl;
62             y.write();
63         }
66         volScalarField yPlus
67         (
68             IOobject
69             (
70                 "yPlus",
71                 runTime.timeName(),
72                 mesh,
73                 IOobject::NO_READ,
74                 IOobject::NO_WRITE
75             ),
76             mesh,
77             dimensionedScalar("yPlus", dimless, 0.0)
78         );
80         Info<< "Reading field U\n" << endl;
81         volVectorField U
82         (
83             IOobject
84             (
85                 "U",
86                 runTime.timeName(),
87                 mesh,
88                 IOobject::MUST_READ,
89                 IOobject::NO_WRITE
90             ),
91             mesh
92         );
94 #       include "createPhi.H"
96         singlePhaseTransportModel laminarTransport(U, phi);
98         autoPtr<incompressible::LESModel> sgsModel
99         (
100             incompressible::LESModel::New(U, phi, laminarTransport)
101         );
103         volScalarField::GeometricBoundaryField d = nearWallDist(mesh).y();
104         volScalarField nuEff = sgsModel->nuEff();
106         const fvPatchList& patches = mesh.boundary();
108         forAll(patches, patchi)
109         {
110             const fvPatch& currPatch = patches[patchi];
112             if (currPatch.isWall())
113             {
114                 yPlus.boundaryField()[patchi] =
115                     d[patchi]
116                    *sqrt
117                     (
118                         nuEff.boundaryField()[patchi]
119                        *mag(U.boundaryField()[patchi].snGrad())
120                     )
121                    /sgsModel->nu().boundaryField()[patchi];
122                 const scalarField& Yp = yPlus.boundaryField()[patchi];
124                 Info<< "Patch " << patchi
125                     << " named " << currPatch.name()
126                     << " y+ : min: " << min(Yp) << " max: " << max(Yp)
127                     << " average: " << average(Yp) << nl << endl;
128             }
129         }
131         Info<< "Writing yPlus to field "
132             << yPlus.name() << nl << endl;
134         yPlus.write();
135     }
137     Info<< "End\n" << endl;
139     return 0;
143 // ************************************************************************* //