Go to the documentation of this file.
40 namespace functionObjects
50 void Foam::functionObjects::yPlus::writeFileHeader(Ostream&
os)
const
74 useWallFunction_(
true)
78 writeFileHeader(file());
87 mesh_.time().timeName(),
97 mesh_.objectRegistry::store(yPlusPtr);
107 useWallFunction_ =
dict.getOrDefault(
"useWallFunction",
true);
118 auto&
yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
122 volScalarField::Boundary& yPlusBf =
yPlus.boundaryFieldRef();
125 lookupObject<turbulenceModel>
131 const volScalarField::Boundary& d = nwd.
y();
135 const volScalarField::Boundary& nutBf = tnut().boundaryField();
148 isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi])
158 yPlusBf[patchi] = nutPf.
yPlus();
160 else if (isA<wallFvPatch>(
patch))
175 <<
"Unable to find turbulence model in the "
176 <<
"database: yPlus will not be calculated" <<
endl;
181 <<
"Please try to use the solver option -postProcess, e.g.:"
182 <<
" <solver> -postProcess -func yPlus" <<
endl;
201 const volScalarField::Boundary& yPlusBf =
yPlus.boundaryField();
208 if (isA<wallFvPatch>(
patch))
212 const scalar minYplus =
gMin(yPlusp);
213 const scalar maxYplus =
gMax(yPlusp);
214 const scalar avgYplus =
gAverage(yPlusp);
219 <<
" y+ : min = " << minYplus <<
", max = " << maxYplus
220 <<
", average = " << avgYplus <<
nl;
222 writeCurrentTime(file());
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
virtual bool read(const dictionary &)
Read the yPlus data.
A class for handling words, derived from Foam::string.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Type gAverage(const FieldField< Field, Type > &f)
bool read(const char *buf, int32_t &val)
Same as readInt32.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
static bool master(const label communicator=worldComm)
Am I the master process.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
#define forAll(list, i)
Loop across all elements in list.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
virtual bool read(const dictionary &dict)
Read.
virtual bool write()
Write the yPlus field.
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool read(const dictionary &dict)
Read optional controls.
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Computes the near wall for turbulence models.
virtual bool execute()
Calculate the yPlus field.
const word & name() const noexcept
Return the name of this functionObject.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
defineTypeNameAndDebug(ObukhovLength, 0)
const polyBoundaryMesh & patches
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
const volScalarField::Boundary & y() const
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Base class for writing single files from the function objects.
Type gMin(const FieldField< Field, Type > &f)
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
#define WarningInFunction
Report a warning using Foam::Warning.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Type gMax(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionSet dimless
Dimensionless.
const volVectorField & U() const
Access function to velocity field.