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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "yPlus.H"
30 #include "volFields.H"
31 #include "turbulenceModel.H"
33 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
43  addToRunTimeSelectionTable(functionObject, yPlus, dictionary);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::functionObjects::yPlus::writeFileHeader(Ostream& os) const
51 {
52  writeHeader(os, "y+ ()");
53 
54  writeCommented(os, "Time");
55  writeTabbed(os, "patch");
56  writeTabbed(os, "min");
57  writeTabbed(os, "max");
58  writeTabbed(os, "average");
59  os << endl;
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const word& name,
68  const Time& runTime,
69  const dictionary& dict
70 )
71 :
73  writeFile(obr_, name, typeName, dict),
74  useWallFunction_(true),
75  writeFields_(true) // May change in the future
76 {
77  read(dict);
78 
79  writeFileHeader(file());
80 
81  volScalarField* yPlusPtr
82  (
83  new volScalarField
84  (
85  IOobject
86  (
87  scopedName(typeName),
88  mesh_.time().timeName(),
89  mesh_,
93  ),
94  mesh_,
96  )
97  );
98 
99  mesh_.objectRegistry::store(yPlusPtr);
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 {
108  {
109  useWallFunction_ = true;
110  writeFields_ = true; // May change in the future
111 
112  dict.readIfPresent("useWallFunction", useWallFunction_);
113  dict.readIfPresent("writeFields", writeFields_);
114 
115  return true;
116  }
117 
118  return false;
119 }
120 
121 
123 {
124  auto& yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
125 
126  if (foundObject<turbulenceModel>(turbulenceModel::propertiesName))
127  {
128  volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
129 
130  const turbulenceModel& model =
131  lookupObject<turbulenceModel>
132  (
134  );
135 
136  const nearWallDist nwd(mesh_);
137  const volScalarField::Boundary& d = nwd.y();
138 
139  // nut needed for wall function patches
140  tmp<volScalarField> tnut = model.nut();
141  const volScalarField::Boundary& nutBf = tnut().boundaryField();
142 
143  // U needed for plain wall patches
144  const volVectorField::Boundary& UBf = model.U().boundaryField();
145 
146  const fvPatchList& patches = mesh_.boundary();
147 
148  forAll(patches, patchi)
149  {
150  const fvPatch& patch = patches[patchi];
151 
152  const auto* nutWallPatch =
153  isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]);
154 
155  if (useWallFunction_ && nutWallPatch)
156  {
157  yPlusBf[patchi] = nutWallPatch->yPlus();
158  }
159  else if (isA<wallFvPatch>(patch))
160  {
161  yPlusBf[patchi] =
162  d[patchi]
163  *sqrt
164  (
165  model.nuEff(patchi)
166  *mag(UBf[patchi].snGrad())
167  )/model.nu(patchi);
168  }
169  }
170  }
171  else
172  {
174  << "Unable to find turbulence model in the "
175  << "database: yPlus will not be calculated" << endl;
176 
177  if (postProcess)
178  {
180  << "Please try to use the solver option -postProcess, e.g.:"
181  << " <solver> -postProcess -func yPlus" << endl;
182  }
183 
184  return false;
185  }
186 
187  return true;
188 }
189 
190 
192 {
193  const auto& yPlus = obr_.lookupObject<volScalarField>(scopedName(typeName));
194 
195  Log << type() << ' ' << name() << " write:" << nl;
196 
197  if (writeFields_)
198  {
199  Log << " writing field " << yPlus.name() << endl;
200  yPlus.write();
201  }
202 
203  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
204  const fvPatchList& patches = mesh_.boundary();
205 
206  forAll(patches, patchi)
207  {
208  const fvPatch& patch = patches[patchi];
209 
210  if (isA<wallFvPatch>(patch))
211  {
212  const scalarField& yPlusp = yPlusBf[patchi];
213 
214  const scalar minYplus = gMin(yPlusp);
215  const scalar maxYplus = gMax(yPlusp);
216  const scalar avgYplus = gAverage(yPlusp);
217 
218  if (UPstream::master())
219  {
220  writeCurrentTime(file());
221  file()
222  << token::TAB << patch.name()
223  << token::TAB << minYplus
224  << token::TAB << maxYplus
225  << token::TAB << avgYplus
226  << endl;
227  }
228 
229  Log << " patch " << patch.name()
230  << " y+ : min = " << minYplus
231  << ", max = " << maxYplus
232  << ", average = " << avgYplus << endl;
233  }
234  }
235 
236  return true;
237 }
238 
239 
240 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:98
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
virtual bool write()
Report min/max/avg (per patch) and log to file, write the yPlus volume field.
Definition: yPlus.C:184
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:115
Type gMin(const FieldField< Field, Type > &f)
Tab [isspace].
Definition: token.H:129
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Abstract base class for turbulence models (RAS, LES and laminar).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest...
Definition: nearWallDist.H:49
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
Type gMax(const FieldField< Field, Type > &f)
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: yPlus.C:59
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
scalar yPlus
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Base class for writing single files from the function objects.
Definition: writeFile.H:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329