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-2018 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 "turbulenceModel.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 turbulenceModel& turbModel
68 )
69 {
70  tmp<volScalarField> tyPlus
71  (
72  new volScalarField
73  (
74  IOobject
75  (
76  type(),
77  mesh_.time().timeName(),
78  mesh_
79  ),
80  mesh_,
81  dimensionedScalar("0", dimless, 0.0)
82  )
83  );
84 
85  volScalarField::Boundary& yPlusBf = tyPlus.ref().boundaryFieldRef();
86 
87  volScalarField::Boundary d = nearWallDist(mesh_).y();
88 
89  const volScalarField::Boundary nutBf =
90  turbModel.nut()().boundaryField();
91 
92  const volScalarField::Boundary nuEffBf =
93  turbModel.nuEff()().boundaryField();
94 
95  const volScalarField::Boundary nuBf =
96  turbModel.nu()().boundaryField();
97 
98  const fvPatchList& patches = mesh_.boundary();
99 
100  forAll(patches, patchi)
101  {
102  const fvPatch& patch = patches[patchi];
103 
104  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
105  {
106  const nutWallFunctionFvPatchScalarField& nutPf =
107  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
108  (
109  nutBf[patchi]
110  );
111 
112  yPlusBf[patchi] = nutPf.yPlus();
113  }
114  else if (isA<wallFvPatch>(patch))
115  {
116  yPlusBf[patchi] =
117  d[patchi]
118  *sqrt
119  (
120  nuEffBf[patchi]
121  *mag(turbModel.U().boundaryField()[patchi].snGrad())
122  )/nuBf[patchi];
123  }
124  }
125 
126  return tyPlus;
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131 
132 Foam::functionObjects::yPlus::yPlus
133 (
134  const word& name,
135  const Time& runTime,
136  const dictionary& dict
137 )
138 :
139  fvMeshFunctionObject(name, runTime, dict),
140  logFiles(obr_, name),
141  writeLocalObjects(obr_, log)
142 {
143  read(dict);
144  resetName(typeName);
145  resetLocalObjectName(typeName);
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
158 {
161 
162  return true;
163 }
164 
165 
167 {
168  if (mesh_.foundObject<turbulenceModel>(turbulenceModel::propertiesName))
169  {
170  const turbulenceModel& model = mesh_.lookupObject<turbulenceModel>
171  (
173  );
174 
175  word name(type());
176 
177  return store(name, calcYPlus(model));
178  }
179  else
180  {
182  << "Unable to find turbulence model in the "
183  << "database" << exit(FatalError);
184  }
185 
186  return true;
187 }
188 
189 
191 {
192  Log << type() << " " << name() << " write:" << nl;
193 
195 
196  logFiles::write();
197 
198  const volScalarField& yPlus =
199  mesh_.lookupObject<volScalarField>(type());
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  writeTime(file());
223  file()
224  << tab << patch.name()
225  << tab << minYplus
226  << tab << maxYplus
227  << tab << avgYplus
228  << endl;
229  }
230  }
231  }
232 
233  Log << endl;
234 
235  return true;
236 }
237 
238 
239 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:157
virtual bool write()
Write function.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:190
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
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:264
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:166
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
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:421
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
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:68
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~yPlus()
Destructor.
Definition: yPlus.C:151
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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:63
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.
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
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;...
Definition: yPlus.H:95
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
A class for managing temporary objects.
Definition: PtrList.H:53
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57