yPlus.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  yPlus
26 
27 Description
28  Calculates and reports yPlus for the near-wall cells of all wall patches,
29  for the specified times for laminar, LES and RAS.
30 
31  For walls at which wall-functions are applied the wall-function provides
32  the y+ values otherwise they are obtained directly from the near-wall
33  velocity gradient and effective and laminar viscosities.
34 
35  Compressible modes is automatically selected based on the existence of the
36  "thermophysicalProperties" dictionary required to construct the
37  thermodynamics package.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
46 #include "nearWallDist.H"
47 #include "wallFvPatch.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 template<class TurbulenceModel>
52 void calcYPlus
53 (
54  const TurbulenceModel& turbulenceModel,
55  const fvMesh& mesh,
56  const volVectorField& U,
58 )
59 {
60  volScalarField::GeometricBoundaryField d = nearWallDist(mesh).y();
61 
62  const volScalarField::GeometricBoundaryField nutBf =
63  turbulenceModel->nut()().boundaryField();
64 
65  const volScalarField::GeometricBoundaryField nuEffBf =
66  turbulenceModel->nuEff()().boundaryField();
67 
68  const volScalarField::GeometricBoundaryField nuBf =
70 
71  const fvPatchList& patches = mesh.boundary();
72 
74  {
75  const fvPatch& patch = patches[patchi];
76 
77  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
78  {
79  const nutWallFunctionFvPatchScalarField& nutPf =
80  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
81  (
82  nutBf[patchi]
83  );
84 
85  yPlus.boundaryField()[patchi] = nutPf.yPlus();
86  const scalarField& Yp = yPlus.boundaryField()[patchi];
87 
88  Info<< "Patch " << patchi
89  << " named " << nutPf.patch().name()
90  << ", wall-function " << nutPf.type()
91  << ", y+ : min: " << gMin(Yp) << " max: " << gMax(Yp)
92  << " average: " << gAverage(Yp) << nl << endl;
93  }
94  else if (isA<wallFvPatch>(patch))
95  {
96  yPlus.boundaryField()[patchi] =
97  d[patchi]
98  *sqrt
99  (
100  nuEffBf[patchi]
101  *mag(U.boundaryField()[patchi].snGrad())
102  )/nuBf[patchi];
103  const scalarField& Yp = yPlus.boundaryField()[patchi];
104 
105  Info<< "Patch " << patchi
106  << " named " << patch.name()
107  << " y+ : min: " << gMin(Yp) << " max: " << gMax(Yp)
108  << " average: " << gAverage(Yp) << nl << endl;
109  }
110  }
111 }
112 
113 
114 void calcIncompressibleYPlus
115 (
116  const fvMesh& mesh,
117  const Time& runTime,
118  const volVectorField& U,
120 )
121 {
122  #include "createPhi.H"
123 
124  singlePhaseTransportModel laminarTransport(U, phi);
125 
126  autoPtr<incompressible::turbulenceModel> turbulenceModel
127  (
129  );
130 
131  calcYPlus(turbulenceModel, mesh, U, yPlus);
132 }
133 
134 
135 void calcCompressibleYPlus
136 (
137  const fvMesh& mesh,
138  const Time& runTime,
139  const volVectorField& U,
141 )
142 {
143  IOobject rhoHeader
144  (
145  "rho",
146  runTime.timeName(),
147  mesh,
148  IOobject::MUST_READ,
149  IOobject::NO_WRITE
150  );
151 
152  if (!rhoHeader.headerOk())
153  {
154  Info<< " no rho field" << endl;
155  return;
156  }
157 
158  Info<< "Reading field rho\n" << endl;
159  volScalarField rho(rhoHeader, mesh);
160 
161  #include "compressibleCreatePhi.H"
162 
163  autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
164  fluidThermo& thermo = pThermo();
165 
166  autoPtr<compressible::turbulenceModel> turbulenceModel
167  (
169  (
170  rho,
171  U,
172  phi,
173  thermo
174  )
175  );
176 
177  calcYPlus(turbulenceModel, mesh, U, yPlus);
178 }
179 
180 
181 int main(int argc, char *argv[])
182 {
183  timeSelector::addOptions();
184  #include "addRegionOption.H"
185  #include "setRootCase.H"
186  #include "createTime.H"
187  instantList timeDirs = timeSelector::select(runTime, args, "yPlus");
188  #include "createNamedMesh.H"
189 
190  forAll(timeDirs, timeI)
191  {
192  runTime.setTime(timeDirs[timeI], timeI);
193  Info<< "Time = " << runTime.timeName() << endl;
194  mesh.readUpdate();
195 
197  (
198  IOobject
199  (
200  "yPlus",
201  runTime.timeName(),
202  mesh,
203  IOobject::NO_READ,
204  IOobject::NO_WRITE
205  ),
206  mesh,
207  dimensionedScalar("yPlus", dimless, 0.0)
208  );
209 
210  IOobject UHeader
211  (
212  "U",
213  runTime.timeName(),
214  mesh,
215  IOobject::MUST_READ,
216  IOobject::NO_WRITE
217  );
218 
219  if (UHeader.headerOk())
220  {
221  Info<< "Reading field U\n" << endl;
223 
224  if
225  (
226  IOobject
227  (
229  runTime.constant(),
230  mesh
231  ).headerOk()
232  )
233  {
234  calcCompressibleYPlus(mesh, runTime, U, yPlus);
235  }
236  else
237  {
238  calcIncompressibleYPlus(mesh, runTime, U, yPlus);
239  }
240  }
241  else
242  {
243  Info<< " no U field" << endl;
244  }
245 
246  Info<< "Writing yPlus to field " << yPlus.name() << nl << endl;
247 
248  yPlus.write();
249  }
250 
251  Info<< "End\n" << endl;
252 
253  return 0;
254 }
255 
256 
257 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
pThermo
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiThermo > pThermo(psiThermo::New(mesh))
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
singlePhaseTransportModel.H
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
turbulentTransportModel.H
nutWallFunctionFvPatchScalarField.H
wallFvPatch.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
U
U
Definition: pEqn.H:46
Foam::fvPatchList
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
dictName
const word dictName("particleTrackDict")
UHeader
IOobject UHeader("U", runTime.timeName(), mesh, IOobject::MUST_READ)
nearWallDist.H
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:60
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
addRegionOption.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
compressibleCreatePhi.H
Creates and initialises the face-flux field phi.
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
scalarField
volScalarField scalarField(fieldObject, mesh)
patchi
label patchi
Definition: getPatchFieldScalar.H:1
createTime.H
fvCFD.H
patches
patches[0]
Definition: createSingleCellMesh.H:36
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
args
Foam::argList args(argc, argv)
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
turbulentFluidThermoModel.H