Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace functionObjects
41{
44}
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void 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 useWallFunction_(true),
75 writeFields_(true) // May change in the future
76{
77 read(dict);
78
79 writeFileHeader(file());
80
81 volScalarField* yPlusPtr
82 (
84 (
86 (
89 mesh_.thisDb(),
93 ),
94 mesh_,
96 )
97 );
99 regIOobject::store(yPlusPtr);
100}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
106{
108 {
109 useWallFunction_ = true;
110 writeFields_ = true; // May change in the future
111
112 dict.readIfPresent("useWallFunction", useWallFunction_);
113 dict.readIfPresent("writeFields", writeFields_);
114
115 return true;
116 }
117
118 return false;
119}
120
121
123{
124 auto& yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
125
126 if (foundObject<turbulenceModel>(turbulenceModel::propertiesName))
127 {
128 volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
129
130 const turbulenceModel& model =
131 lookupObject<turbulenceModel>
132 (
134 );
135
136 const nearWallDist nwd(mesh_);
137 const volScalarField::Boundary& d = nwd.y();
138
139 // nut needed for wall function patches
140 tmp<volScalarField> tnut = model.nut();
141 const volScalarField::Boundary& nutBf = tnut().boundaryField();
142
143 // U needed for plain wall patches
144 const volVectorField::Boundary& UBf = model.U().boundaryField();
145
146 const fvPatchList& patches = mesh_.boundary();
147
148 forAll(patches, patchi)
149 {
150 const fvPatch& patch = patches[patchi];
151
152 const auto* nutWallPatch =
154
155 if (useWallFunction_ && nutWallPatch)
156 {
157 yPlusBf[patchi] = nutWallPatch->yPlus();
158 }
159 else if (isA<wallFvPatch>(patch))
160 {
161 yPlusBf[patchi] =
162 d[patchi]
163 *sqrt
164 (
165 model.nuEff(patchi)
166 *mag(UBf[patchi].snGrad())
167 )/model.nu(patchi);
168 }
169 }
170 }
171 else
172 {
174 << "Unable to find turbulence model in the "
175 << "database: yPlus will not be calculated" << endl;
176
177 if (postProcess)
178 {
180 << "Please try to use the solver option -postProcess, e.g.:"
181 << " <solver> -postProcess -func yPlus" << endl;
182 }
183
184 return false;
185 }
186
187 return true;
188}
189
190
192{
193 const auto& yPlus = obr_.lookupObject<volScalarField>(scopedName(typeName));
194
195 Log << type() << ' ' << name() << " write:" << nl;
196
197 if (writeFields_)
198 {
199 Log << " writing field " << yPlus.name() << endl;
200 yPlus.write();
201 }
202
203 const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
204 const fvPatchList& patches = mesh_.boundary();
205
206 forAll(patches, patchi)
207 {
208 const fvPatch& patch = patches[patchi];
209
210 if (isA<wallFvPatch>(patch))
211 {
212 const scalarField& yPlusp = yPlusBf[patchi];
213
214 auto limits = gMinMax(yPlusp);
215 auto avg = gAverage(yPlusp);
216
217 if (UPstream::master())
218 {
219 writeCurrentTime(file());
220 file()
221 << token::TAB << patch.name()
222 << token::TAB << limits.min()
223 << token::TAB << limits.max()
224 << token::TAB << avg
225 << endl;
226 }
227
228 Log << " patch " << patch.name()
229 << " y+ : min = " << limits.min()
230 << ", max = " << limits.max()
231 << ", average = " << avg << endl;
232 }
233 }
234
235 return true;
236}
237
238
239// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
static bool postProcess
Global post-processing mode switch.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the magnitude of an input field.
Definition mag.H:136
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
const objectRegistry & obr_
Reference to the region objectRegistry.
ObjectType & lookupObjectRef(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.
Definition writeFile.H:113
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition writeFile.C:344
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
Computes the near wall for turbulence models.
Definition yPlus.H:158
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition yPlus.C:98
yPlus(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition yPlus.C:59
virtual bool execute()
Execute the function-object operations.
Definition yPlus.C:115
virtual bool write()
Write the function-object results.
Definition yPlus.C:184
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
const volScalarField::Boundary & y() const
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
@ TAB
Tab [isspace].
Definition token.H:142
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.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
auto limits
Definition setRDeltaT.H:186
const polyBoundaryMesh & patches
engineTime & runTime
scalar yPlus
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
const dimensionSet dimless
Dimensionless.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299