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  (
43  functionObject,
44  yPlus,
45  dictionary
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 
96  forAll(patches, patchi)
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_.lookupObject<momentumTransportModel>
173  (
175  (
176  momentumTransportModel::typeName,
177  phaseName_
178  )
179  );
180 
181  word name(IOobject::groupName(type(), phaseName_));
182 
183  return store(name, calcYPlus(model));
184  }
185  else
186  {
188  << "Unable to find turbulence model in the "
189  << "database" << exit(FatalError);
190  }
191 
192  return true;
193 }
194 
195 
197 {
198  Log << type() << " " << name() << " write:" << nl;
199 
201 
202  logFiles::write();
203 
204  const volScalarField& yPlus =
206  (
207  IOobject::groupName(type(), phaseName_)
208  );
209 
210  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
211  const fvPatchList& patches = mesh_.boundary();
212 
213  forAll(patches, patchi)
214  {
215  const fvPatch& patch = patches[patchi];
216 
217  if (isA<wallFvPatch>(patch))
218  {
219  const scalarField& yPlusp = yPlusBf[patchi];
220 
221  const scalar minYplus = gMin(yPlusp);
222  const scalar maxYplus = gMax(yPlusp);
223  const scalar avgYplus = gAverage(yPlusp);
224 
225  if (Pstream::master())
226  {
227  Log << " patch " << patch.name()
228  << " y+ : min = " << minYplus << ", max = " << maxYplus
229  << ", average = " << avgYplus << nl;
230 
231  writeTime(file());
232  file()
233  << tab << patch.name()
234  << tab << minYplus
235  << tab << maxYplus
236  << tab << avgYplus
237  << endl;
238  }
239  }
240  }
241 
242  Log << endl;
243 
244  return true;
245 }
246 
247 
248 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:152
virtual bool write()
Write function.
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:103
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:196
virtual bool write()
Write function.
Definition: logFiles.C:167
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:259
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:164
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
const dimensionSet dimless
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
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.
compressibleMomentumTransportModel momentumTransportModel
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: yPlus.C:129
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~yPlus()
Destructor.
Definition: yPlus.C:146
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
scalar yPlus
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Evaluates and outputs turbulence y+ for models. Values written to time directories as field &#39;yPlus&#39; o...
Definition: yPlus.H:101
A class for managing temporary objects.
Definition: PtrList.H:53
const volScalarField::Boundary & y() const
Access to the near-wall distances.
Definition: nearWallDist.H:93
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57