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-2020 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"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
39 
41  (
42  functionObject,
43  yPlus,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::functionObjects::yPlus::writeFileHeader(const label i)
53 {
54  writeHeader(file(), "y+ ()");
55 
56  writeCommented(file(), "Time");
57  writeTabbed(file(), "patch");
58  writeTabbed(file(), "min");
59  writeTabbed(file(), "max");
60  writeTabbed(file(), "average");
61  file() << endl;
62 }
63 
64 
65 Foam::tmp<Foam::volScalarField> Foam::functionObjects::yPlus::calcYPlus
66 (
67  const momentumTransportModel& turbModel
68 )
69 {
70  tmp<volScalarField> tyPlus
71  (
73  (
74  IOobject::groupName(type(), phaseName_),
75  mesh_,
77  )
78  );
79 
80  volScalarField::Boundary& yPlusBf = tyPlus.ref().boundaryFieldRef();
81 
82  volScalarField::Boundary d = nearWallDist(mesh_).y();
83 
84  const volScalarField::Boundary nutBf =
85  turbModel.nut()().boundaryField();
86 
87  const volScalarField::Boundary nuEffBf =
88  turbModel.nuEff()().boundaryField();
89 
90  const volScalarField::Boundary nuBf =
91  turbModel.nu()().boundaryField();
92 
93  const fvPatchList& patches = mesh_.boundary();
94 
95  forAll(patches, patchi)
96  {
97  const fvPatch& patch = patches[patchi];
98 
99  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
100  {
101  const nutWallFunctionFvPatchScalarField& nutPf =
102  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
103  (
104  nutBf[patchi]
105  );
106 
107  yPlusBf[patchi] = nutPf.yPlus();
108  }
109  else if (isA<wallFvPatch>(patch))
110  {
111  yPlusBf[patchi] =
112  d[patchi]
113  *sqrt
114  (
115  nuEffBf[patchi]
116  *mag(turbModel.U().boundaryField()[patchi].snGrad())
117  )/nuBf[patchi];
118  }
119  }
120 
121  return tyPlus;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 (
129  const word& name,
130  const Time& runTime,
131  const dictionary& dict
132 )
133 :
134  fvMeshFunctionObject(name, runTime, dict),
135  logFiles(obr_, name),
136  writeLocalObjects(obr_, log),
137  phaseName_(dict.lookupOrDefault<word>("phase", word::null))
138 {
139  read(dict);
140  resetName(IOobject::groupName(typeName, phaseName_));
141  resetLocalObjectName(IOobject::groupName(typeName, phaseName_));
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146 
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 {
157 
158  return true;
159 }
160 
161 
163 {
164  if (mesh_.foundObject<momentumTransportModel>
165  (
166  IOobject::groupName(momentumTransportModel::typeName, phaseName_))
167  )
168  {
169  const momentumTransportModel& model =
170  mesh_.lookupObject<momentumTransportModel>
171  (
173  (
174  momentumTransportModel::typeName,
175  phaseName_
176  )
177  );
178 
179  word name(IOobject::groupName(type(), phaseName_));
180 
181  return store(name, calcYPlus(model));
182  }
183  else
184  {
186  << "Unable to find turbulence model in the "
187  << "database" << exit(FatalError);
188  }
189 
190  return true;
191 }
192 
193 
195 {
196  Log << type() << " " << name() << " write:" << nl;
197 
199 
200  logFiles::write();
201 
202  const volScalarField& yPlus =
204  (
205  IOobject::groupName(type(), phaseName_)
206  );
207 
208  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
209  const fvPatchList& patches = mesh_.boundary();
210 
211  forAll(patches, patchi)
212  {
213  const fvPatch& patch = patches[patchi];
214 
215  if (isA<wallFvPatch>(patch))
216  {
217  const scalarField& yPlusp = yPlusBf[patchi];
218 
219  const scalar minYplus = gMin(yPlusp);
220  const scalar maxYplus = gMax(yPlusp);
221  const scalar avgYplus = gAverage(yPlusp);
222 
223  if (Pstream::master())
224  {
225  Log << " patch " << patch.name()
226  << " y+ : min = " << minYplus << ", max = " << maxYplus
227  << ", average = " << avgYplus << nl;
228 
229  writeTime(file());
230  file()
231  << tab << patch.name()
232  << tab << minYplus
233  << tab << maxYplus
234  << tab << avgYplus
235  << endl;
236  }
237  }
238  }
239 
240  Log << endl;
241 
242  return true;
243 }
244 
245 
246 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:153
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:194
virtual bool write()
Write function.
Definition: logFiles.C:180
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
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:162
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:61
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
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:68
Macros for easy insertion into run-time selection tables.
CompressibleMomentumTransportModel< fluidThermo > 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:128
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:147
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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:143
Specialization 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
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57