yPlus.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "yPlus.H"
27 #include "momentumTransportModel.H"
29 #include "wallFvPatch.H"
30 #include "nearWallDist.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
40 
42  (
44  yPlus,
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::functionObjects::yPlus::writeFileHeader(const label i)
54 {
55  writeHeader(file(), "y+ ()");
56 
57  writeCommented(file(), "Time");
58  writeTabbed(file(), "patch");
59  writeTabbed(file(), "min");
60  writeTabbed(file(), "max");
61  writeTabbed(file(), "average");
62  file() << endl;
63 }
64 
65 
66 Foam::tmp<Foam::volScalarField> Foam::functionObjects::yPlus::calcYPlus
67 (
68  const momentumTransportModel& turbModel
69 )
70 {
71  tmp<volScalarField> tyPlus
72  (
74  (
75  IOobject::groupName(type(), phaseName_),
76  mesh_,
78  )
79  );
80 
81  volScalarField::Boundary& yPlusBf = tyPlus.ref().boundaryFieldRef();
82 
83  const volScalarField::Boundary& d = nearWallDist::New(mesh_).y();
84 
85  const tmp<volScalarField> tnut = turbModel.nut();
86  const volScalarField::Boundary& nutBf = tnut().boundaryField();
87 
88  const tmp<volScalarField> tnuEff = turbModel.nuEff();
89  const volScalarField::Boundary& nuEffBf = tnuEff().boundaryField();
90 
91  const tmp<volScalarField> tnu = turbModel.nu();
92  const volScalarField::Boundary& nuBf = tnu().boundaryField();
93 
94  const fvPatchList& patches = mesh_.boundary();
95 
97  {
98  const fvPatch& patch = patches[patchi];
99 
100  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
101  {
102  const nutWallFunctionFvPatchScalarField& nutPf =
103  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
104  (
105  nutBf[patchi]
106  );
107 
108  yPlusBf[patchi] = nutPf.yPlus();
109  }
110  else if (isA<wallFvPatch>(patch))
111  {
112  yPlusBf[patchi] =
113  d[patchi]
114  *sqrt
115  (
116  nuEffBf[patchi]
117  *mag(turbModel.U().boundaryField()[patchi].snGrad())
118  )/nuBf[patchi];
119  }
120  }
121 
122  return tyPlus;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const word& name,
131  const Time& runTime,
132  const dictionary& dict
133 )
134 :
135  fvMeshFunctionObject(name, runTime, dict),
136  logFiles(obr_, name),
137  writeLocalObjects(obr_, log),
138  phaseName_(dict.lookupOrDefault<word>("phase", word::null))
139 {
140  read(dict);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
156 
157  resetName(IOobject::groupName(typeName, phaseName_));
158  resetLocalObjectName(IOobject::groupName(typeName, phaseName_));
159 
160  return true;
161 }
162 
163 
165 {
166  if (mesh_.foundObject<momentumTransportModel>
167  (
168  IOobject::groupName(momentumTransportModel::typeName, phaseName_))
169  )
170  {
171  const momentumTransportModel& model =
172  mesh_.lookupType<momentumTransportModel>(phaseName_);
173 
174  word name(IOobject::groupName(type(), phaseName_));
175 
176  return store(name, calcYPlus(model));
177  }
178  else
179  {
181  << "Unable to find turbulence model in the "
182  << "database" << exit(FatalError);
183  }
184 
185  return true;
186 }
187 
188 
190 {
191  Log << type() << " " << name() << " write:" << nl;
192 
194 
195  logFiles::write();
196 
197  const volScalarField& yPlus =
199  (
200  IOobject::groupName(type(), phaseName_)
201  );
202 
203  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
204  const fvPatchList& patches = mesh_.boundary();
205 
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 (Pstream::master())
219  {
220  Log << " patch " << patch.name()
221  << " y+ : min = " << minYplus << ", max = " << maxYplus
222  << ", average = " << avgYplus << nl;
223 
224  writeTime(file());
225  file()
226  << tab << patch.name()
227  << tab << minYplus
228  << tab << maxYplus
229  << tab << avgYplus
230  << endl;
231  }
232  }
233  }
234 
235  Log << endl;
236 
237  return true;
238 }
239 
240 
241 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static nearWallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricBoundaryField class.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:60
OFstream & file()
Return access to the file (if only 1)
Definition: logFiles.C:113
virtual bool write()
Write function.
Definition: logFiles.C:173
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:106
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
virtual bool read(const dictionary &)
Read optional controls.
void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:121
void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:131
void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:110
FunctionObject base class for managing a list of objects on behalf of the inheriting function object,...
virtual bool read(const dictionary &)
Read the list of objects to be written.
virtual bool write()
Write function.
Evaluates and outputs turbulence y+ for models. Values written to time directories as field 'yPlus' o...
Definition: yPlus.H:106
virtual ~yPlus()
Destructor.
Definition: yPlus.C:146
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: yPlus.C:129
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:164
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:189
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:152
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & y() const
Access to the near-wall distances.
Definition: nearWallDist.H:99
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const fvPatchList & patches
#define Log
Report write to Foam::Info if the local log switch is true.
compressibleMomentumTransportModel momentumTransportModel
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char tab
Definition: Ostream.H:265
const dimensionSet dimless
dimensionedScalar sqrt(const dimensionedScalar &ds)
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
Type gMin(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:266
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict