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-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
39namespace Foam
40{
41namespace functionObjects
42{
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 : patchSet_)
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:
92 writeFile(mesh_, name, typeName, dict),
93 patchSet_()
94{
95 read(dict);
96
98
99 volVectorField* wallShearStressPtr
100 (
102 (
104 (
105 scopedName(typeName),
106 mesh_.time().timeName(),
107 mesh_,
110 ),
111 mesh_,
113 )
114 );
115
116 mesh_.objectRegistry::store(wallShearStressPtr);
117}
118
119
120// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121
123{
126
127 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
128
129 patchSet_ =
130 mesh_.boundaryMesh().patchSet
131 (
132 dict.getOrDefault<wordRes>("patches", wordRes())
133 );
134
135 Info<< type() << " " << name() << ":" << nl;
136
137 if (patchSet_.empty())
138 {
139 forAll(pbm, patchi)
140 {
141 if (isA<wallPolyPatch>(pbm[patchi]))
142 {
143 patchSet_.insert(patchi);
144 }
145 }
146
147 Info<< " processing all wall patches" << nl << endl;
148 }
149 else
150 {
151 Info<< " processing wall patches: " << nl;
152 labelHashSet filteredPatchSet;
153 for (const label patchi : patchSet_)
154 {
155 if (isA<wallPolyPatch>(pbm[patchi]))
156 {
157 filteredPatchSet.insert(patchi);
158 Info<< " " << pbm[patchi].name() << endl;
159 }
160 else
161 {
163 << "Requested wall shear stress on non-wall boundary "
164 << "type patch: " << pbm[patchi].name() << endl;
165 }
166 }
167
168 Info<< endl;
169
170 patchSet_ = filteredPatchSet;
171 }
172
173 return true;
174}
175
176
178{
179 auto& wallShearStress =
180 mesh_.lookupObjectRef<volVectorField>(scopedName(typeName));
181
182 // Compressible
183 {
184 typedef compressible::turbulenceModel turbType;
185
186 const turbType* modelPtr =
187 findObject<turbType>(turbulenceModel::propertiesName);
188
189 if (modelPtr)
190 {
191 calcShearStress(modelPtr->devRhoReff(), wallShearStress);
192 return true;
193 }
194 }
195
196 // Incompressible
197 {
198 typedef incompressible::turbulenceModel turbType;
199
200 const turbType* modelPtr =
201 findObject<turbType>(turbulenceModel::propertiesName);
202
203 if (modelPtr)
204 {
205 calcShearStress(modelPtr->devReff(), wallShearStress);
206 return true;
207 }
208 }
209
211 << "Unable to find turbulence model in the "
212 << "database" << exit(FatalError);
213
214 return false;
215}
216
217
219{
220 const auto& wallShearStress =
221 obr_.lookupObject<volVectorField>(scopedName(typeName));
222
223 Log << type() << " " << name() << " write:" << nl
224 << " writing field " << wallShearStress.name() << endl;
225
227
228 const fvPatchList& patches = mesh_.boundary();
229
230 for (const label patchi : patchSet_)
231 {
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 writeCurrentTime(file());
242
243 file()
244 << token::TAB << pp.name()
245 << token::TAB << minSsp
246 << token::TAB << maxSsp
247 << endl;
248 }
249
250 Log << " min/max(" << pp.name() << ") = "
251 << minSsp << ", " << maxSsp << endl;
252 }
253
254 return true;
255}
256
257
258// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const dimensionSet & dimensions() const
Return dimensions.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
Templated abstract base class for single-phase incompressible turbulence models.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:149
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
ObjectType & lookupObjectRef(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Computes the wall-shear stress at selected wall patches.
virtual void writeFileHeader(Ostream &os) const
File header information.
virtual bool execute()
Calculate the wall shear-stress.
void calcShearStress(const volSymmTensorField &Reff, volVectorField &shearStress)
Calculate the shear-stress.
virtual bool write()
Write the wall shear-stress.
virtual bool read(const dictionary &)
Read the wallShearStress data.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:285
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:233
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
splitCell * master() const
Definition: splitCell.H:113
static const word propertiesName
Default name of the turbulence properties dictionary.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.