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-2023 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 volSymmTensorField& 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  label patchi = iter.key();
82 
83  const vectorField& Sfp = mesh_.Sf().boundaryField()[patchi];
84  const scalarField& magSfp = mesh_.magSf().boundaryField()[patchi];
85 
86  wallShearStressBf[patchi] = (-Sfp/magSfp) & tau.boundaryField()[patchi];
87  }
88 
89  return twallShearStress;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
98  const Time& runTime,
99  const dictionary& dict
100 )
101 :
102  fvMeshFunctionObject(name, runTime, dict),
103  logFiles(obr_, name),
104  writeLocalObjects(obr_, log),
105  phaseName_(word::null),
106  patchSet_()
107 {
108  read(dict);
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
114 
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
125 
126  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
127 
128  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
129 
130  patchSet_ =
131  mesh_.boundaryMesh().patchSet
132  (
133  wordReList(dict.lookupOrDefault("patches", wordReList()))
134  );
135 
136  Info<< type() << " " << name() << ":" << nl;
137 
138  if (patchSet_.empty())
139  {
140  forAll(pbm, patchi)
141  {
142  if (isA<wallPolyPatch>(pbm[patchi]))
143  {
144  patchSet_.insert(patchi);
145  }
146  }
147 
148  Info<< " processing all wall patches" << nl << endl;
149  }
150  else
151  {
152  Info<< " processing wall patches: " << nl;
153  labelHashSet filteredPatchSet;
154  forAllConstIter(labelHashSet, patchSet_, iter)
155  {
156  label patchi = iter.key();
157  if (isA<wallPolyPatch>(pbm[patchi]))
158  {
159  filteredPatchSet.insert(patchi);
160  Info<< " " << pbm[patchi].name() << endl;
161  }
162  else
163  {
165  << "Requested wall shear stress on non-wall boundary "
166  << "type patch: " << pbm[patchi].name() << endl;
167  }
168  }
169 
170  Info<< endl;
171 
172  patchSet_ = filteredPatchSet;
173  }
174 
175  resetName(typeName);
176  resetLocalObjectName(typeName);
177 
178  return true;
179 }
180 
181 
183 {
184  const word fieldName(IOobject::groupName(type(), phaseName_));
185 
186  typedef compressible::momentumTransportModel cmpModel;
188 
189  const word momentumTransportModelName
190  (
191  IOobject::groupName(momentumTransportModel::typeName, phaseName_)
192  );
193 
194  if (mesh_.foundObject<cmpModel>(momentumTransportModelName))
195  {
196  const cmpModel& model =
197  mesh_.lookupObject<cmpModel>(momentumTransportModelName);
198 
199  return store(fieldName, calcShearStress(model.devTau()));
200  }
201  else if (mesh_.foundObject<icoModel>(momentumTransportModelName))
202  {
203  const icoModel& model =
204  mesh_.lookupObject<icoModel>(momentumTransportModelName);
205 
206  return store(fieldName, calcShearStress(model.devSigma()));
207  }
208  else
209  {
211  << "Unable to find turbulence model in the "
212  << "database" << exit(FatalError);
213 
214  return false;
215  }
216 }
217 
218 
220 {
221  Log << type() << " " << name() << " write:" << nl;
222 
224 
225  logFiles::write();
226 
228  (
229  IOobject::groupName(type(), phaseName_)
230  );
231 
232  const fvPatchList& patches = mesh_.boundary();
233 
234  forAllConstIter(labelHashSet, patchSet_, iter)
235  {
236  label patchi = iter.key();
237  const fvPatch& pp = patches[patchi];
238 
239  const vectorField& ssp = wallShearStress.boundaryField()[patchi];
240 
241  vector minSsp = gMin(ssp);
242  vector maxSsp = gMax(ssp);
243 
244  if (Pstream::master())
245  {
246  file() << time_.userTimeValue()
247  << tab << pp.name()
248  << tab << minSsp
249  << tab << maxSsp
250  << endl;
251  }
252 
253  Log << " min/max(" << pp.name() << ") = "
254  << minSsp << ", " << maxSsp << endl;
255  }
256 
257  Log << endl;
258 
259  return true;
260 }
261 
262 
263 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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:107
virtual bool write()
Write function.
Definition: logFiles.C:167
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:145
Base class for single-phase incompressible turbulence models.
Foam::polyBoundaryMesh.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
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:306
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
#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:251
static const char tab
Definition: Ostream.H:259
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
error FatalError
Type gMin(const FieldField< Field, Type > &f)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:64
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
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.