wallShearStress.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallShearStress.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
31 #include "wallPolyPatch.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
47 
48 void Foam::functionObjects::wallShearStress::writeFileHeader(const label i)
49 {
50  // Add headers to output data
51  writeHeader(file(), "Wall shear stress");
52  writeCommented(file(), "Time");
53  writeTabbed(file(), "patch");
54  writeTabbed(file(), "min");
55  writeTabbed(file(), "max");
56  file() << endl;
57 }
58 
59 
61 Foam::functionObjects::wallShearStress::calcShearStress
62 (
63  const surfaceVectorField& tau
64 )
65 {
66  tmp<volVectorField> twallShearStress
67  (
69  (
70  type(),
71  mesh_,
72  dimensionedVector(tau.dimensions(), Zero)
73  )
74  );
75 
76  volVectorField::Boundary& wallShearStressBf =
77  twallShearStress.ref().boundaryFieldRef();
78 
79  forAllConstIter(labelHashSet, patchSet_, iter)
80  {
81  const label patchi = iter.key();
82  wallShearStressBf[patchi] = -tau.boundaryField()[patchi];
83  }
84 
85  return twallShearStress;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const word& name,
94  const Time& runTime,
95  const dictionary& dict
96 )
97 :
98  fvMeshFunctionObject(name, runTime, dict),
99  logFiles(obr_, name),
100  writeLocalObjects(obr_, log),
101  phaseName_(word::null),
102  patchSet_()
103 {
104  read(dict);
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110 
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
121 
122  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
123 
124  patchSet_ = mesh_.boundaryMesh().patchSet(dict, true);
125 
126  Info<< type() << " " << name() << ":" << nl;
127 
128  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
129 
130  if (patchSet_.empty())
131  {
132  forAll(pbm, patchi)
133  {
134  if (isA<wallPolyPatch>(pbm[patchi]))
135  {
136  patchSet_.insert(patchi);
137  }
138  }
139 
140  Info<< " processing all wall patches" << nl << endl;
141  }
142  else
143  {
144  Info<< " processing wall patches: " << nl;
145 
146  labelHashSet filteredPatchSet;
147  forAllConstIter(labelHashSet, patchSet_, iter)
148  {
149  const label patchi = iter.key();
150  filteredPatchSet.insert(patchi);
151 
152  if (isA<wallPolyPatch>(pbm[patchi]))
153  {
154  Info<< " " << pbm[patchi].name() << endl;
155  }
156  else
157  {
158  Info<< " " << pbm[patchi].name()
159  << " type " << pbm[patchi].type() << endl;
160  }
161  }
162 
163  Info<< endl;
164 
165  patchSet_ = filteredPatchSet;
166  }
167 
168  resetName(typeName);
169  resetLocalObjectName(typeName);
170 
171  return true;
172 }
173 
174 
176 {
177  const word fieldName(IOobject::groupName(type(), phaseName_));
178 
179  typedef compressible::momentumTransportModel cmpModel;
181 
182  const word momentumTransportModelName
183  (
184  IOobject::groupName(momentumTransportModel::typeName, phaseName_)
185  );
186 
187  if (mesh_.foundObject<cmpModel>(momentumTransportModelName))
188  {
189  const cmpModel& model =
190  mesh_.lookupObject<cmpModel>(momentumTransportModelName);
191 
192  store(fieldName, calcShearStress(model.devTau()));
193 
194  return true;
195  }
196  else if (mesh_.foundObject<icoModel>(momentumTransportModelName))
197  {
198  const icoModel& model =
199  mesh_.lookupObject<icoModel>(momentumTransportModelName);
200 
201  store(fieldName, calcShearStress(model.devSigma()));
202 
203  return true;
204  }
205  else
206  {
208  << "Unable to find turbulence model in the "
209  << "database" << exit(FatalError);
210 
211  return false;
212  }
213 }
214 
215 
217 {
218  Log << type() << " " << name() << " write:" << nl;
219 
221 
222  logFiles::write();
223 
225  (
226  IOobject::groupName(type(), phaseName_)
227  );
228 
229  const fvPatchList& patches = mesh_.boundary();
230 
231  forAllConstIter(labelHashSet, patchSet_, iter)
232  {
233  label patchi = iter.key();
234  const fvPatch& pp = patches[patchi];
235 
236  const vectorField& ssp = wallShearStress.boundaryField()[patchi];
237 
238  vector minSsp = gMin(ssp);
239  vector maxSsp = gMax(ssp);
240 
241  if (Pstream::master())
242  {
243  file() << time_.userTimeValue()
244  << tab << pp.name()
245  << tab << minSsp
246  << tab << maxSsp
247  << endl;
248  }
249 
250  Log << " min/max(" << pp.name() << ") = "
251  << minSsp << ", " << maxSsp << endl;
252  }
253 
254  Log << endl;
255 
256  return true;
257 }
258 
259 
260 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:60
OFstream & file()
Return access to the file (if only 1)
Definition: logFiles.C:113
virtual bool write()
Write function.
Definition: logFiles.C:173
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:106
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
virtual bool read(const dictionary &)
Read optional controls.
Calculates and write the shear-stress at wall patches as the volVectorField field 'wallShearStress' o...
wallShearStress(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
virtual bool execute()
Calculate the wall shear-stress.
virtual bool write()
Write the wall shear-stress.
virtual bool read(const dictionary &)
Read the wallShearStress data.
void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:121
void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:131
void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:110
FunctionObject base class for managing a list of objects on behalf of the inheriting function object,...
void resetLocalObjectName(const word &name)
Reset the list of local object names from a single word.
virtual bool read(const dictionary &)
Read the list of objects to be written.
virtual bool write()
Write function.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
Base class for single-phase incompressible turbulence models.
Foam::polyBoundaryMesh.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const vector tau
label patchi
const fvPatchList & patches
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
static const char tab
Definition: Ostream.H:266
messageStream Info
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
error FatalError
Type gMin(const FieldField< Field, Type > &f)
SurfaceField< vector > surfaceVectorField
static const char nl
Definition: Ostream.H:267
Type gMax(const FieldField< Field, Type > &f)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
Foam::surfaceFields.