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 | Copyright (C) 2013-2016 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 
61 (
62  const fvMesh& mesh,
63  const volSymmTensorField& Reff,
64  volVectorField& shearStress
65 )
66 {
67  forAllConstIter(labelHashSet, patchSet_, iter)
68  {
69  label patchi = iter.key();
70 
71  vectorField& ssp = shearStress.boundaryFieldRef()[patchi];
72  const vectorField& Sfp = mesh.Sf().boundaryField()[patchi];
73  const scalarField& magSfp = mesh.magSf().boundaryField()[patchi];
74  const symmTensorField& Reffp = Reff.boundaryField()[patchi];
75 
76  ssp = (-Sfp/magSfp) & Reffp;
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
83 Foam::functionObjects::wallShearStress::wallShearStress
84 (
85  const word& name,
86  const Time& runTime,
87  const dictionary& dict
88 )
89 :
90  writeFiles(name, runTime, dict, name),
91  patchSet_()
92 {
93  if (!isA<fvMesh>(obr_))
94  {
96  << "objectRegistry is not an fvMesh" << exit(FatalError);
97  }
98 
99  const fvMesh& mesh = refCast<const fvMesh>(obr_);
100 
101  volVectorField* wallShearStressPtr
102  (
103  new volVectorField
104  (
105  IOobject
106  (
107  type(),
108  mesh.time().timeName(),
109  mesh,
112  ),
113  mesh,
115  (
116  "0",
118  Zero
119  )
120  )
121  );
122 
123  mesh.objectRegistry::store(wallShearStressPtr);
124 
125  read(dict);
126  resetName(typeName);
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
131 
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  writeFiles::read(dict);
141 
142  const fvMesh& mesh = refCast<const fvMesh>(obr_);
143  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
144 
145  patchSet_ =
146  mesh.boundaryMesh().patchSet
147  (
148  wordReList(dict.lookupOrDefault("patches", wordReList()))
149  );
150 
151  Info<< type() << " " << name() << ":" << nl;
152 
153  if (patchSet_.empty())
154  {
155  forAll(pbm, patchi)
156  {
157  if (isA<wallPolyPatch>(pbm[patchi]))
158  {
159  patchSet_.insert(patchi);
160  }
161  }
162 
163  Info<< " processing all wall patches" << nl << endl;
164  }
165  else
166  {
167  Info<< " processing wall patches: " << nl;
168  labelHashSet filteredPatchSet;
169  forAllConstIter(labelHashSet, patchSet_, iter)
170  {
171  label patchi = iter.key();
172  if (isA<wallPolyPatch>(pbm[patchi]))
173  {
174  filteredPatchSet.insert(patchi);
175  Info<< " " << pbm[patchi].name() << endl;
176  }
177  else
178  {
180  << "Requested wall shear stress on non-wall boundary "
181  << "type patch: " << pbm[patchi].name() << endl;
182  }
183  }
184 
185  Info<< endl;
186 
187  patchSet_ = filteredPatchSet;
188  }
189 
190  return true;
191 }
192 
193 
195 {
196  typedef compressible::turbulenceModel cmpModel;
197  typedef incompressible::turbulenceModel icoModel;
198 
199  const fvMesh& mesh = refCast<const fvMesh>(obr_);
200 
202  const_cast<volVectorField&>
203  (
205  );
206 
208  if (mesh.foundObject<cmpModel>(turbulenceModel::propertiesName))
209  {
210  const cmpModel& model =
212 
213  Reff = model.devRhoReff();
214  }
215  else if (mesh.foundObject<icoModel>(turbulenceModel::propertiesName))
216  {
217  const icoModel& model =
219 
220  Reff = model.devReff();
221  }
222  else
223  {
225  << "Unable to find turbulence model in the "
226  << "database" << exit(FatalError);
227  }
228 
229  calcShearStress(mesh, Reff(), wallShearStress);
230 
231  return true;
232 }
233 
234 
236 {
238 
241 
242  Log << type() << " " << name() << " write:" << nl
243  << " writing field " << wallShearStress.name() << endl;
244 
245  wallShearStress.write();
246 
247  const fvMesh& mesh = refCast<const fvMesh>(obr_);
248  const fvPatchList& patches = mesh.boundary();
249 
250  forAllConstIter(labelHashSet, patchSet_, iter)
251  {
252  label patchi = iter.key();
253  const fvPatch& pp = patches[patchi];
254 
255  const vectorField& ssp = wallShearStress.boundaryField()[patchi];
256 
257  vector minSsp = gMin(ssp);
258  vector maxSsp = gMax(ssp);
259 
260  if (Pstream::master())
261  {
262  file() << mesh.time().value()
263  << token::TAB << pp.name()
264  << token::TAB << minSsp
265  << token::TAB << maxSsp
266  << endl;
267  }
268 
269  Log << " min/max(" << pp.name() << ") = "
270  << minSsp << ", " << maxSsp << endl;
271  }
272 
273  return true;
274 }
275 
276 
277 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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 surfaceVectorField & Sf() const
Return cell face area vectors.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual bool read(const dictionary &)
Read the wallShearStress data.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool foundObject(const word &name) const
Is the named Type found?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
void calcShearStress(const fvMesh &mesh, const volSymmTensorField &Reff, volVectorField &shearStress)
Calculate the shear stress.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
functionObject base class for writing files
Definition: writeFiles.H:56
const word & name() const
Return name.
Definition: fvPatch.H:149
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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:411
const Type & value() const
Return const reference to value.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
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
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
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 Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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.
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
This function object evaluates and outputs the shear stress at wall patches. The result is written as...
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.
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:262
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual bool write()
Write the wall shear-stress.
label patchi
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual bool write()
Write function.
Definition: writeFiles.C:197
virtual bool write() const
Write using setting from DB.
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
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:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
const ObjectType & lookupObject(const word &fieldName) const
Lookup field from the objectRegistry.
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243