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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
38namespace Foam
39{
40namespace functionObjects
41{
44}
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void 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{
76 read(dict);
77
78 writeFileHeader(file());
79
80 volScalarField* yPlusPtr
81 (
83 (
85 (
86 scopedName(typeName),
88 mesh_,
91 ),
92 mesh_,
94 )
95 );
96
97 mesh_.objectRegistry::store(yPlusPtr);
98}
99
100
101// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102
104{
106 {
107 useWallFunction_ = dict.getOrDefault("useWallFunction", true);
108
109 return true;
110 }
111
112 return false;
113}
114
115
117{
118 auto& yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
119
120 if (foundObject<turbulenceModel>(turbulenceModel::propertiesName))
121 {
122 volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
123
124 const turbulenceModel& model =
125 lookupObject<turbulenceModel>
126 (
128 );
129
130 const nearWallDist nwd(mesh_);
131 const volScalarField::Boundary& d = nwd.y();
132
133 // nut needed for wall function patches
134 tmp<volScalarField> tnut = model.nut();
135 const volScalarField::Boundary& nutBf = tnut().boundaryField();
136
137 // U needed for plain wall patches
138 const volVectorField::Boundary& UBf = model.U().boundaryField();
139
140 const fvPatchList& patches = mesh_.boundary();
141
142 forAll(patches, patchi)
143 {
144 const fvPatch& patch = patches[patchi];
145
146 if
147 (
148 isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi])
149 && useWallFunction_
150 )
151 {
153 dynamic_cast<const nutWallFunctionFvPatchScalarField&>
154 (
155 nutBf[patchi]
156 );
157
158 yPlusBf[patchi] = nutPf.yPlus();
159 }
160 else if (isA<wallFvPatch>(patch))
161 {
162 yPlusBf[patchi] =
163 d[patchi]
164 *sqrt
165 (
166 model.nuEff(patchi)
167 *mag(UBf[patchi].snGrad())
168 )/model.nu(patchi);
169 }
170 }
171 }
172 else
173 {
175 << "Unable to find turbulence model in the "
176 << "database: yPlus will not be calculated" << endl;
177
178 if (postProcess)
179 {
181 << "Please try to use the solver option -postProcess, e.g.:"
182 << " <solver> -postProcess -func yPlus" << endl;
183 }
184
185 return false;
186 }
187
188 return true;
189}
190
191
193{
194 const auto& yPlus = obr_.lookupObject<volScalarField>(scopedName(typeName));
195
196 Log << type() << " " << name() << " write:" << nl
197 << " writing field " << yPlus.name() << endl;
198
199 yPlus.write();
200
201 const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
202 const fvPatchList& patches = mesh_.boundary();
203
204 forAll(patches, patchi)
205 {
206 const fvPatch& patch = patches[patchi];
207
208 if (isA<wallFvPatch>(patch))
209 {
210 const scalarField& yPlusp = yPlusBf[patchi];
211
212 const scalar minYplus = gMin(yPlusp);
213 const scalar maxYplus = gMax(yPlusp);
214 const scalar avgYplus = gAverage(yPlusp);
215
216 if (Pstream::master())
217 {
218 Log << " patch " << patch.name()
219 << " y+ : min = " << minYplus << ", max = " << maxYplus
220 << ", average = " << avgYplus << nl;
221
222 writeCurrentTime(file());
223 file()
224 << token::TAB << patch.name()
225 << token::TAB << minYplus
226 << token::TAB << maxYplus
227 << token::TAB << avgYplus
228 << endl;
229 }
230 }
231 }
232
233 return true;
234}
235
236
237// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
Definition: mag.H:153
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:285
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:233
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
Computes the near wall for turbulence models.
Definition: yPlus.H:160
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:116
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:192
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:103
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:56
const volScalarField::Boundary & y() const
Definition: nearWallDist.H:89
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333