40namespace functionObjects
50void Foam::functionObjects::yPlus::writeFileHeader(Ostream&
os)
const
74 useWallFunction_(true)
78 writeFileHeader(
file());
97 mesh_.objectRegistry::store(yPlusPtr);
107 useWallFunction_ =
dict.getOrDefault(
"useWallFunction",
true);
118 auto&
yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
125 lookupObject<turbulenceModel>
148 isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi])
158 yPlusBf[patchi] = nutPf.
yPlus();
160 else if (isA<wallFvPatch>(patch))
167 *
mag(UBf[patchi].snGrad())
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;
208 if (isA<wallFvPatch>(patch))
212 const scalar minYplus =
gMin(yPlusp);
213 const scalar maxYplus =
gMax(yPlusp);
214 const scalar avgYplus =
gAverage(yPlusp);
218 Log <<
" patch " << patch.name()
219 <<
" y+ : min = " << minYplus <<
", max = " << maxYplus
220 <<
", average = " << avgYplus <<
nl;
222 writeCurrentTime(file());
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static word timeName(const scalar t, const int precision=precision_)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Base class for writing single files from the function objects.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
virtual OFstream & file()
Return access to the file (if only 1)
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Computes the near wall for turbulence models.
virtual bool execute()
Calculate the yPlus field.
virtual bool write()
Write the yPlus field.
virtual bool read(const dictionary &)
Read the yPlus data.
const Time & time() const
Return the top-level database.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
const volScalarField::Boundary & y() const
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
splitCell * master() const
A class for managing temporary objects.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const polyBoundaryMesh & patches
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimless
Dimensionless.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.