yPlus.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2013-2017 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 {
38  defineTypeNameAndDebug(yPlus, 0);
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 void Foam::functionObjects::yPlus::calcYPlus
66 (
67  const turbulenceModel& turbModel,
68  volScalarField& yPlus
69 )
70 {
71  volScalarField::Boundary d = nearWallDist(mesh_).y();
72 
73  const volScalarField::Boundary nutBf =
74  turbModel.nut()().boundaryField();
75 
76  const volScalarField::Boundary nuEffBf =
77  turbModel.nuEff()().boundaryField();
78 
79  const volScalarField::Boundary nuBf =
80  turbModel.nu()().boundaryField();
81 
82  const fvPatchList& patches = mesh_.boundary();
83 
84  volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
85 
86  forAll(patches, patchi)
87  {
88  const fvPatch& patch = patches[patchi];
89 
90  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
91  {
92  const nutWallFunctionFvPatchScalarField& nutPf =
93  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
94  (
95  nutBf[patchi]
96  );
97 
98  yPlusBf[patchi] = nutPf.yPlus();
99  }
100  else if (isA<wallFvPatch>(patch))
101  {
102  yPlusBf[patchi] =
103  d[patchi]
104  *sqrt
105  (
106  nuEffBf[patchi]
107  *mag(turbModel.U().boundaryField()[patchi].snGrad())
108  )/nuBf[patchi];
109  }
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
115 
116 Foam::functionObjects::yPlus::yPlus
117 (
118  const word& name,
119  const Time& runTime,
120  const dictionary& dict
121 )
122 :
123  fvMeshFunctionObject(name, runTime, dict),
124  logFiles(obr_, name),
125  writeLocalObjects(obr_, log)
126 {
127  volScalarField* yPlusPtr
128  (
129  new volScalarField
130  (
131  IOobject
132  (
133  type(),
134  mesh_.time().timeName(),
135  mesh_,
138  ),
139  mesh_,
140  dimensionedScalar("0", dimless, 0.0)
141  )
142  );
143 
144  mesh_.objectRegistry::store(yPlusPtr);
145 
146  read(dict);
147  resetName(typeName);
148  resetLocalObjectName(typeName);
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
164 
165  return true;
166 }
167 
168 
170 {
171  volScalarField& yPlus =
172  mesh_.lookupObjectRef<volScalarField>(type());
173 
174  if (mesh_.foundObject<turbulenceModel>(turbulenceModel::propertiesName))
175  {
176  const turbulenceModel& model = mesh_.lookupObject<turbulenceModel>
177  (
179  );
180 
181  calcYPlus(model, yPlus);
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 =
203  mesh_.lookupObject<volScalarField>(type());
204 
205  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
206  const fvPatchList& patches = mesh_.boundary();
207 
208  forAll(patches, patchi)
209  {
210  const fvPatch& patch = patches[patchi];
211 
212  if (isA<wallFvPatch>(patch))
213  {
214  const scalarField& yPlusp = yPlusBf[patchi];
215 
216  const scalar minYplus = gMin(yPlusp);
217  const scalar maxYplus = gMax(yPlusp);
218  const scalar avgYplus = gAverage(yPlusp);
219 
220  if (Pstream::master())
221  {
222  Log << " patch " << patch.name()
223  << " y+ : min = " << minYplus << ", max = " << maxYplus
224  << ", average = " << avgYplus << nl;
225 
226  writeTime(file());
227  file()
228  << token::TAB << patch.name()
229  << token::TAB << minYplus
230  << token::TAB << maxYplus
231  << token::TAB << avgYplus
232  << endl;
233  }
234  }
235  }
236 
237  Log << endl;
238 
239  return true;
240 }
241 
242 
243 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:160
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: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
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:169
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:253
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:412
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:262
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:485
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~yPlus()
Destructor.
Definition: yPlus.C:154
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
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...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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