wallShearStress.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) 2017-2023 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 "wallShearStress.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
34 #include "wallPolyPatch.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(wallShearStress, 0);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
50 
52 {
53  writeHeader(os, "Wall shear stress");
54  writeCommented(os, "Time");
55  writeTabbed(os, "patch");
56  writeTabbed(os, "min");
57  writeTabbed(os, "max");
58  os << endl;
59 }
60 
61 
63 (
64  const volSymmTensorField& Reff,
65  volVectorField& shearStress
66 )
67 {
68  shearStress.dimensions().reset(Reff.dimensions());
69 
70  for (const label patchi : patchIDs_)
71  {
72  vectorField& ssp = shearStress.boundaryFieldRef()[patchi];
73  const vectorField& Sfp = mesh_.Sf().boundaryField()[patchi];
74  const scalarField& magSfp = mesh_.magSf().boundaryField()[patchi];
75  const symmTensorField& Reffp = Reff.boundaryField()[patchi];
76 
77  ssp = (-Sfp/magSfp) & Reffp;
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 (
86  const word& name,
87  const Time& runTime,
88  const dictionary& dict
89 )
90 :
91  fvMeshFunctionObject(name, runTime, dict),
92  writeFile(mesh_, name, typeName, dict),
93  writeFields_(true) // May change in the future
94 {
95  read(dict);
96 
98 
99  volVectorField* wallShearStressPtr
100  (
101  new volVectorField
102  (
103  IOobject
104  (
105  scopedName(typeName),
106  mesh_.time().timeName(),
107  mesh_,
111  ),
112  mesh_,
114  )
115  );
117  mesh_.objectRegistry::store(wallShearStressPtr);
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 {
127 
128  writeFields_ = true; // May change in the future
129  dict.readIfPresent("writeFields", writeFields_);
130 
131  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
132 
134  labelHashSet patchSet(0);
135  if (dict.readIfPresent("patches", patchNames) && !patchNames.empty())
136  {
137  patchSet = pbm.patchSet(patchNames);
138  }
139 
140  labelHashSet allWalls(pbm.findPatchIDs<wallPolyPatch>());
141 
142  Info<< type() << ' ' << name() << ':' << nl;
143 
144  if (patchSet.empty())
145  {
146  patchIDs_ = allWalls.sortedToc();
147 
148  Info<< " processing all (" << patchIDs_.size()
149  << ") wall patches" << nl << endl;
150  }
151  else
152  {
153  allWalls &= patchSet;
154  patchSet -= allWalls;
155  patchIDs_ = allWalls.sortedToc();
156 
157  if (!patchSet.empty())
158  {
160  << "Requested wall shear stress on ("
161  << patchSet.size() << ") non-wall patches:" << nl;
162 
163  for (const label patchi : patchSet.sortedToc())
164  {
165  Info<< " " << pbm[patchi].name() << nl;
166  }
167  Info<< nl;
168  }
169 
170  Info<< " processing (" << patchIDs_.size()
171  << ") wall patches:" << nl;
172 
173  for (const label patchi : patchIDs_)
174  {
175  Info<< " " << pbm[patchi].name() << nl;
176  }
178  }
179 
180  return true;
181 }
182 
183 
185 {
186  auto& wallShearStress =
187  mesh_.lookupObjectRef<volVectorField>(scopedName(typeName));
188 
189  // Compressible
190  {
191  typedef compressible::turbulenceModel turbType;
192 
193  const turbType* modelPtr =
194  findObject<turbType>(turbulenceModel::propertiesName);
195 
196  if (modelPtr)
197  {
198  calcShearStress(modelPtr->devRhoReff(), wallShearStress);
199  return true;
200  }
201  }
202 
203  // Incompressible
204  {
205  typedef incompressible::turbulenceModel turbType;
206 
207  const turbType* modelPtr =
208  findObject<turbType>(turbulenceModel::propertiesName);
209 
210  if (modelPtr)
211  {
212  calcShearStress(modelPtr->devReff(), wallShearStress);
213  return true;
214  }
215  }
216 
218  << "Unable to find turbulence model in the "
219  << "database" << exit(FatalError);
220 
221  return false;
222 }
223 
224 
226 {
227  const auto& wallShearStress =
228  obr_.lookupObject<volVectorField>(scopedName(typeName));
229 
230  Log << type() << ' ' << name() << " write:" << nl;
231 
232  if (writeFields_)
233  {
234  Log << " writing field " << wallShearStress.name() << endl;
236  }
237 
238  const fvPatchList& patches = mesh_.boundary();
239 
240  for (const label patchi : patchIDs_)
241  {
242  const fvPatch& pp = patches[patchi];
243 
244  const vectorField& ssp = wallShearStress.boundaryField()[patchi];
245 
246  const MinMax<vector> limits = gMinMax(ssp);
247 
248  if (UPstream::master())
249  {
250  writeCurrentTime(file());
251 
252  file()
253  << token::TAB << pp.name()
254  << token::TAB << limits.min()
255  << token::TAB << limits.max()
256  << endl;
257  }
258 
259  Log << " min/max(" << pp.name() << ") = "
260  << limits.min() << ", " << limits.max() << endl;
261  }
262 
263  return true;
264 }
265 
266 
267 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & pbm
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual bool read(const dictionary &)
Read the wallShearStress data.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Tab [isspace].
Definition: token.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A min/max value pair with additional methods. In addition to conveniently storing values...
Definition: HashSet.H:72
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
wallShearStress(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Ignore writing from objectRegistry::writeObject()
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Abstract base-class for Time/database function objects.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool execute()
Calculate the wall shear-stress.
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
void calcShearStress(const volSymmTensorField &Reff, volVectorField &shearStress)
Calculate the shear-stress.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:46
wordList patchNames(nPatches)
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
Computes the wall-shear stress at selected wall patches.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
virtual void writeFileHeader(Ostream &os) const
File header information.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:142
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:337
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
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
virtual bool write()
Report min/max and log to file, write the wall shear-stress volume field.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
#define Log
Definition: PDRblock.C:28
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:121
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329