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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "yPlus.H"
30 #include "volFields.H"
31 #include "turbulenceModel.H"
33 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
43  addToRunTimeSelectionTable(functionObject, yPlus, dictionary);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::functionObjects::yPlus::writeFileHeader(Ostream& os) const
51 {
52  writeHeader(os, "y+ ()");
53 
54  writeCommented(os, "Time");
55  writeTabbed(os, "patch");
56  writeTabbed(os, "min");
57  writeTabbed(os, "max");
58  writeTabbed(os, "average");
59  os << endl;
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const word& name,
68  const Time& runTime,
69  const dictionary& dict
70 )
71 :
73  writeFile(obr_, name, typeName, dict)
74 {
75  read(dict);
76 
77  writeFileHeader(file());
78 
79  volScalarField* yPlusPtr
80  (
81  new volScalarField
82  (
83  IOobject
84  (
85  typeName,
86  mesh_.time().timeName(),
87  mesh_,
90  ),
91  mesh_,
93  )
94  );
95 
96  mesh_.objectRegistry::store(yPlusPtr);
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 {
106 
107  return true;
108 }
109 
110 
112 {
114  lookupObjectRef<volScalarField>(typeName);
115 
116  if (foundObject<turbulenceModel>(turbulenceModel::propertiesName))
117  {
118  volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
119 
120  const turbulenceModel& model =
121  lookupObject<turbulenceModel>
122  (
124  );
125 
126  const nearWallDist nwd(mesh_);
127  const volScalarField::Boundary& d = nwd.y();
128 
129  // nut needed for wall function patches
130  tmp<volScalarField> tnut = model.nut();
131  const volScalarField::Boundary& nutBf = tnut().boundaryField();
132 
133  // U needed for plain wall patches
134  const volVectorField::Boundary& UBf = model.U().boundaryField();
135 
136  const fvPatchList& patches = mesh_.boundary();
137 
138  forAll(patches, patchi)
139  {
140  const fvPatch& patch = patches[patchi];
141 
142  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
143  {
144  const nutWallFunctionFvPatchScalarField& nutPf =
145  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
146  (
147  nutBf[patchi]
148  );
149 
150  yPlusBf[patchi] = nutPf.yPlus();
151  }
152  else if (isA<wallFvPatch>(patch))
153  {
154  yPlusBf[patchi] =
155  d[patchi]
156  *sqrt
157  (
158  model.nuEff(patchi)
159  *mag(UBf[patchi].snGrad())
160  )/model.nu(patchi);
161  }
162  }
163  }
164  else
165  {
167  << "Unable to find turbulence model in the "
168  << "database: yPlus will not be calculated" << endl;
169  return false;
170  }
171 
172  return true;
173 }
174 
175 
177 {
178  const volScalarField& yPlus =
180 
181  Log << type() << " " << name() << " write:" << nl
182  << " writing field " << yPlus.name() << endl;
183 
184  yPlus.write();
185 
186  const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
187  const fvPatchList& patches = mesh_.boundary();
188 
189  forAll(patches, patchi)
190  {
191  const fvPatch& patch = patches[patchi];
192 
193  if (isA<wallFvPatch>(patch))
194  {
195  const scalarField& yPlusp = yPlusBf[patchi];
196 
197  const scalar minYplus = gMin(yPlusp);
198  const scalar maxYplus = gMax(yPlusp);
199  const scalar avgYplus = gAverage(yPlusp);
200 
201  if (Pstream::master())
202  {
203  Log << " patch " << patch.name()
204  << " y+ : min = " << minYplus << ", max = " << maxYplus
205  << ", average = " << avgYplus << nl;
206 
207  writeCurrentTime(file());
208  file()
209  << token::TAB << patch.name()
210  << token::TAB << minYplus
211  << token::TAB << maxYplus
212  << token::TAB << avgYplus
213  << endl;
214  }
215  }
216  }
217 
218  return true;
219 }
220 
221 
222 // ************************************************************************* //
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Log
#define Log
Definition: PDRblock.C:35
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::functionObjects::yPlus::read
virtual bool read(const dictionary &)
Read the yPlus data.
Definition: yPlus.C:102
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::regionFunctionObject::lookupObject
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Definition: regionFunctionObjectTemplates.C:87
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
nutWallFunctionFvPatchScalarField.H
Foam::nutWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
wallFvPatch.H
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::functionObjects::writeFile::writeHeader
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:299
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:214
yPlus.H
Foam::Field< scalar >
Foam::functionObjects::yPlus::write
virtual bool write()
Write the yPlus field.
Definition: yPlus.C:176
Foam::functionObjects::yPlus::yPlus
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: yPlus.C:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::functionObjects::yPlus
Computes the near wall for turbulence models.
Definition: yPlus.H:137
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::functionObjects::yPlus::execute
virtual bool execute()
Calculate the yPlus field.
Definition: yPlus.C:111
Foam::functionObjects::writeFile::writeCommented
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:273
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::token::TAB
Tab [isspace].
Definition: token.H:118
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::nearWallDist
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:53
Foam::nearWallDist::y
const volScalarField::Boundary & y() const
Definition: nearWallDist.H:89
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::turbulenceModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::IOobject::NO_READ
Definition: IOobject.H:123
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::functionObjects::writeFile::writeTabbed
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:289
turbulenceModel.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144