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  (
73  (
74  type(),
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 {
138  read(dict);
139  resetName(typeName);
140  resetLocalObjectName(typeName);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
156 
157  return true;
158 }
159 
160 
162 {
163  if (mesh_.foundObject<turbulenceModel>(turbulenceModel::propertiesName))
164  {
165  const turbulenceModel& model = mesh_.lookupObject<turbulenceModel>
166  (
168  );
169 
170  word name(type());
171 
172  return store(name, calcYPlus(model));
173  }
174  else
175  {
177  << "Unable to find turbulence model in the "
178  << "database" << exit(FatalError);
179  }
180 
181  return true;
182 }
183 
184 
186 {
187  Log << type() << " " << name() << " write:" << nl;
188 
190 
191  logFiles::write();
192 
193  const volScalarField& yPlus =
194  mesh_.lookupObject<volScalarField>(type());
195 
196  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
197  const fvPatchList& patches = mesh_.boundary();
198 
199  forAll(patches, patchi)
200  {
201  const fvPatch& patch = patches[patchi];
202 
203  if (isA<wallFvPatch>(patch))
204  {
205  const scalarField& yPlusp = yPlusBf[patchi];
206 
207  const scalar minYplus = gMin(yPlusp);
208  const scalar maxYplus = gMax(yPlusp);
209  const scalar avgYplus = gAverage(yPlusp);
210 
211  if (Pstream::master())
212  {
213  Log << " patch " << patch.name()
214  << " y+ : min = " << minYplus << ", max = " << maxYplus
215  << ", average = " << avgYplus << nl;
216 
217  writeTime(file());
218  file()
219  << tab << patch.name()
220  << tab << minYplus
221  << tab << maxYplus
222  << tab << avgYplus
223  << endl;
224  }
225  }
226  }
227 
228  Log << endl;
229 
230  return true;
231 }
232 
233 
234 // ************************************************************************* //
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:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:185
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:264
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:161
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
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:423
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
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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)
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)
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~yPlus()
Destructor.
Definition: yPlus.C:146
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: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;...
Definition: yPlus.H:95
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