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-2018 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 {
40  defineTypeNameAndDebug(wallShearStress, 0);
41  addToRunTimeSelectionTable(functionObject, wallShearStress, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
47 
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 
62 (
63  const volSymmTensorField& Reff
64 )
65 {
66  tmp<volVectorField> twallShearStress
67  (
68  new volVectorField
69  (
70  IOobject
71  (
72  type(),
73  mesh_.time().timeName(),
74  mesh_
75  ),
76  mesh_,
77  dimensionedVector("0", Reff.dimensions(), Zero)
78  )
79  );
80 
81  volVectorField::Boundary& wallShearStressBf =
82  twallShearStress.ref().boundaryFieldRef();
83 
84  forAllConstIter(labelHashSet, patchSet_, iter)
85  {
86  label patchi = iter.key();
87 
88  const vectorField& Sfp = mesh_.Sf().boundaryField()[patchi];
89  const scalarField& magSfp = mesh_.magSf().boundaryField()[patchi];
90  const symmTensorField& Reffp = Reff.boundaryField()[patchi];
91 
92  wallShearStressBf[patchi] = (-Sfp/magSfp) & Reffp;
93  }
94 
95  return twallShearStress;
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
101 Foam::functionObjects::wallShearStress::wallShearStress
102 (
103  const word& name,
104  const Time& runTime,
105  const dictionary& dict
106 )
107 :
108  fvMeshFunctionObject(name, runTime, dict),
109  logFiles(obr_, name),
110  writeLocalObjects(obr_, log),
111  patchSet_()
112 {
113  read(dict);
114  resetName(typeName);
115  resetLocalObjectName(typeName);
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
120 
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
128 {
131 
132  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
133 
134  patchSet_ =
135  mesh_.boundaryMesh().patchSet
136  (
137  wordReList(dict.lookupOrDefault("patches", wordReList()))
138  );
139 
140  Info<< type() << " " << name() << ":" << nl;
141 
142  if (patchSet_.empty())
143  {
144  forAll(pbm, patchi)
145  {
146  if (isA<wallPolyPatch>(pbm[patchi]))
147  {
148  patchSet_.insert(patchi);
149  }
150  }
151 
152  Info<< " processing all wall patches" << nl << endl;
153  }
154  else
155  {
156  Info<< " processing wall patches: " << nl;
157  labelHashSet filteredPatchSet;
158  forAllConstIter(labelHashSet, patchSet_, iter)
159  {
160  label patchi = iter.key();
161  if (isA<wallPolyPatch>(pbm[patchi]))
162  {
163  filteredPatchSet.insert(patchi);
164  Info<< " " << pbm[patchi].name() << endl;
165  }
166  else
167  {
169  << "Requested wall shear stress on non-wall boundary "
170  << "type patch: " << pbm[patchi].name() << endl;
171  }
172  }
173 
174  Info<< endl;
175 
176  patchSet_ = filteredPatchSet;
177  }
178 
179  return true;
180 }
181 
182 
184 {
185  typedef compressible::turbulenceModel cmpModel;
186  typedef incompressible::turbulenceModel icoModel;
187 
189  if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
190  {
191  const cmpModel& model =
192  mesh_.lookupObject<cmpModel>(turbulenceModel::propertiesName);
193 
194  Reff = model.devRhoReff();
195  }
196  else if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
197  {
198  const icoModel& model =
199  mesh_.lookupObject<icoModel>(turbulenceModel::propertiesName);
200 
201  Reff = model.devReff();
202  }
203  else
204  {
206  << "Unable to find turbulence model in the "
207  << "database" << exit(FatalError);
208  }
209 
210  word name(type());
211 
212  return store(name, calcShearStress(Reff));
213 }
214 
215 
217 {
218  Log << type() << " " << name() << " write:" << nl;
219 
221 
222  logFiles::write();
223 
226 
227  const fvPatchList& patches = mesh_.boundary();
228 
229  forAllConstIter(labelHashSet, patchSet_, iter)
230  {
231  label patchi = iter.key();
232  const fvPatch& pp = patches[patchi];
233 
234  const vectorField& ssp = wallShearStress.boundaryField()[patchi];
235 
236  vector minSsp = gMin(ssp);
237  vector maxSsp = gMax(ssp);
238 
239  if (Pstream::master())
240  {
241  file() << mesh_.time().value()
242  << tab << pp.name()
243  << tab << minSsp
244  << tab << maxSsp
245  << endl;
246  }
247 
248  Log << " min/max(" << pp.name() << ") = "
249  << minSsp << ", " << maxSsp << endl;
250  }
251 
252  Log << endl;
253 
254  return true;
255 }
256 
257 
258 // ************************************************************************* //
Foam::surfaceFields.
virtual bool write()
Write function.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool write()
Write function.
Definition: logFiles.C:180
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
const word & name() const
Return name.
Definition: IOobject.H:297
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual bool read(const dictionary &)
Read the wallShearStress data.
static const char tab
Definition: Ostream.H:264
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:421
patches[0]
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
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.
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
virtual bool execute()
Calculate the wall shear-stress.
Templated abstract base class for single-phase incompressible turbulence models.
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual void writeFileHeader(const label i)
File header information.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
Calculates and write the shear-stress at wall patches as the volVectorField field &#39;wallShearStress&#39;...
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual bool write()
Write the wall shear-stress.
tmp< volVectorField > calcShearStress(const volSymmTensorField &Reff)
Calculate the shear-stress.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
FunctionObject base class for managing a list of objects on behalf of the inheriting function object...
virtual bool read(const dictionary &)
Read the list of objects to be written.
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57