streamFunction.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) 2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "streamFunction.H"
30 #include "surfaceFields.H"
31 #include "pointFields.H"
32 #include "emptyPolyPatch.H"
33 #include "symmetryPlanePolyPatch.H"
34 #include "symmetryPolyPatch.H"
35 #include "wedgePolyPatch.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace functionObjects
43 {
44  defineTypeNameAndDebug(streamFunction, 0);
45  addToRunTimeSelectionTable(functionObject, streamFunction, dictionary);
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 Foam::tmp<Foam::pointScalarField> Foam::functionObjects::streamFunction::calc
53 (
54  const surfaceScalarField& phi
55 ) const
56 {
57  Log << " functionObjects::" << type() << " " << name()
58  << " calculating stream-function" << endl;
59 
60  Vector<label> slabNormal((Vector<label>::one - mesh_.geometricD())/2);
61  const direction slabDir
62  (
63  slabNormal
65  );
66 
67  const pointMesh& pMesh = pointMesh::New(mesh_);
68 
69  auto tstreamFunction = tmp<pointScalarField>::New
70  (
71  IOobject
72  (
73  "streamFunction",
74  time_.timeName(),
75  mesh_
76  ),
77  pMesh,
78  dimensionedScalar(phi.dimensions(), Zero)
79  );
80  pointScalarField& streamFunction = tstreamFunction.ref();
81 
82  labelList visitedPoint(mesh_.nPoints(), Zero);
83 
84  label nVisited = 0;
85  label nVisitedOld = 0;
86 
87  const faceUList& faces = mesh_.faces();
88  const pointField& points = mesh_.points();
89 
90  label nInternalFaces = mesh_.nInternalFaces();
91 
92  vectorField unitAreas(mesh_.faceAreas());
93  unitAreas /= mag(unitAreas);
94 
96 
97  bool finished = true;
98 
99  // Find the boundary face with zero flux. Set the stream function
100  // to zero on that face
101  bool found = false;
102 
103  do
104  {
105  found = false;
106 
107  forAll(patches, patchi)
108  {
109  const primitivePatch& bouFaces = patches[patchi];
110 
111  if (!isType<emptyPolyPatch>(patches[patchi]))
112  {
113  forAll(bouFaces, facei)
114  {
115  if (magSqr(phi.boundaryField()[patchi][facei]) < SMALL)
116  {
117  const labelList& zeroPoints = bouFaces[facei];
118 
119  // Zero flux face found
120  found = true;
121 
122  forAll(zeroPoints, pointi)
123  {
124  if (visitedPoint[zeroPoints[pointi]] == 1)
125  {
126  found = false;
127  break;
128  }
129  }
130 
131  if (found)
132  {
133  Log << " Zero face: patch: " << patchi
134  << " face: " << facei << endl;
135 
136  forAll(zeroPoints, pointi)
137  {
138  streamFunction[zeroPoints[pointi]] = 0;
139  visitedPoint[zeroPoints[pointi]] = 1;
140  nVisited++;
141  }
142 
143  break;
144  }
145  }
146  }
147  }
148 
149  if (found) break;
150  }
151 
152  if (!found)
153  {
154  Log << " Zero flux boundary face not found. "
155  << "Using cell as a reference."
156  << endl;
157 
158  const cellList& c = mesh_.cells();
159 
160  forAll(c, ci)
161  {
162  labelList zeroPoints = c[ci].labels(mesh_.faces());
163 
164  bool found = true;
165 
166  forAll(zeroPoints, pointi)
167  {
168  if (visitedPoint[zeroPoints[pointi]] == 1)
169  {
170  found = false;
171  break;
172  }
173  }
174 
175  if (found)
176  {
177  forAll(zeroPoints, pointi)
178  {
179  streamFunction[zeroPoints[pointi]] = 0.0;
180  visitedPoint[zeroPoints[pointi]] = 1;
181  nVisited++;
182  }
183 
184  break;
185  }
186  else
187  {
189  << "Cannot find initialisation face or a cell."
190  << exit(FatalError);
191  }
192  }
193  }
194 
195  // Loop through all faces. If one of the points on
196  // the face has the streamfunction value different
197  // from -1, all points with -1 ont that face have the
198  // streamfunction value equal to the face flux in
199  // that point plus the value in the visited point
200  do
201  {
202  finished = true;
203 
204  for (label facei = nInternalFaces; facei<faces.size(); facei++)
205  {
206  const labelList& curBPoints = faces[facei];
207  bool bPointFound = false;
208 
209  scalar currentBStream = 0.0;
210  vector currentBStreamPoint(0, 0, 0);
211 
212  forAll(curBPoints, pointi)
213  {
214  // Check if the point has been visited
215  if (visitedPoint[curBPoints[pointi]] == 1)
216  {
217  // The point has been visited
218  currentBStream = streamFunction[curBPoints[pointi]];
219  currentBStreamPoint = points[curBPoints[pointi]];
220 
221  bPointFound = true;
222 
223  break;
224  }
225  }
226 
227  if (bPointFound)
228  {
229  // Sort out other points on the face
230  forAll(curBPoints, pointi)
231  {
232  // Check if the point has been visited
233  if (visitedPoint[curBPoints[pointi]] == 0)
234  {
235  label patchNo =
236  mesh_.boundaryMesh().whichPatch(facei);
237 
238  if
239  (
240  !isType<emptyPolyPatch>(patches[patchNo])
241  && !isType<symmetryPlanePolyPatch>
242  (patches[patchNo])
243  && !isType<symmetryPolyPatch>(patches[patchNo])
244  && !isType<wedgePolyPatch>(patches[patchNo])
245  )
246  {
247  label faceNo =
248  mesh_.boundaryMesh()[patchNo]
249  .whichFace(facei);
250 
251  vector edgeHat =
252  points[curBPoints[pointi]]
253  - currentBStreamPoint;
254  edgeHat.replace(slabDir, 0);
255  edgeHat.normalise();
256 
257  vector nHat = unitAreas[facei];
258 
259  if (edgeHat.y() > VSMALL)
260  {
261  visitedPoint[curBPoints[pointi]] = 1;
262  nVisited++;
263 
264  streamFunction[curBPoints[pointi]] =
265  currentBStream
266  + phi.boundaryField()[patchNo][faceNo]
267  *sign(nHat.x());
268  }
269  else if (edgeHat.y() < -VSMALL)
270  {
271  visitedPoint[curBPoints[pointi]] = 1;
272  nVisited++;
273 
274  streamFunction[curBPoints[pointi]] =
275  currentBStream
276  - phi.boundaryField()[patchNo][faceNo]
277  *sign(nHat.x());
278  }
279  else
280  {
281  if (edgeHat.x() > VSMALL)
282  {
283  visitedPoint[curBPoints[pointi]] = 1;
284  nVisited++;
285 
286  streamFunction[curBPoints[pointi]] =
287  currentBStream
288  + phi.boundaryField()[patchNo][faceNo]
289  *sign(nHat.y());
290  }
291  else if (edgeHat.x() < -VSMALL)
292  {
293  visitedPoint[curBPoints[pointi]] = 1;
294  nVisited++;
295 
296  streamFunction[curBPoints[pointi]] =
297  currentBStream
298  - phi.boundaryField()[patchNo][faceNo]
299  *sign(nHat.y());
300  }
301  }
302  }
303  }
304  }
305  }
306  else
307  {
308  finished = false;
309  }
310  }
311 
312  for (label facei=0; facei<nInternalFaces; facei++)
313  {
314  // Get the list of point labels for the face
315  const labelList& curPoints = faces[facei];
316 
317  bool pointFound = false;
318 
319  scalar currentStream = 0.0;
320  point currentStreamPoint(0, 0, 0);
321 
322  forAll(curPoints, pointi)
323  {
324  // Check if the point has been visited
325  if (visitedPoint[curPoints[pointi]] == 1)
326  {
327  // The point has been visited
328  currentStream = streamFunction[curPoints[pointi]];
329  currentStreamPoint = points[curPoints[pointi]];
330  pointFound = true;
331 
332  break;
333  }
334  }
335 
336  if (pointFound)
337  {
338  // Sort out other points on the face
339  forAll(curPoints, pointi)
340  {
341  // Check if the point has been visited
342  if (visitedPoint[curPoints[pointi]] == 0)
343  {
344  vector edgeHat =
345  points[curPoints[pointi]] - currentStreamPoint;
346 
347  edgeHat.replace(slabDir, 0);
348  edgeHat.normalise();
349 
350  vector nHat = unitAreas[facei];
351 
352  if (edgeHat.y() > VSMALL)
353  {
354  visitedPoint[curPoints[pointi]] = 1;
355  nVisited++;
356 
357  streamFunction[curPoints[pointi]] =
358  currentStream
359  + phi[facei]*sign(nHat.x());
360  }
361  else if (edgeHat.y() < -VSMALL)
362  {
363  visitedPoint[curPoints[pointi]] = 1;
364  nVisited++;
365 
366  streamFunction[curPoints[pointi]] =
367  currentStream
368  - phi[facei]*sign(nHat.x());
369  }
370  }
371  }
372  }
373  else
374  {
375  finished = false;
376  }
377  }
378 
379  if (nVisited == nVisitedOld)
380  {
381  // Find new seed. This must be a
382  // multiply connected domain
383  Log << " Exhausted a seed, looking for new seed "
384  << "(this is correct for multiply connected domains).";
385 
386  break;
387  }
388  else
389  {
390  nVisitedOld = nVisited;
391  }
392  } while (!finished);
393  } while (!finished);
394 
395  // Normalise the stream-function by the 2D mesh thickness
396  // calculate thickness here to avoid compiler oddness (#2603)
397  const scalar thickness = vector(slabNormal) & mesh_.bounds().span();
398 
399  streamFunction /= thickness;
400  streamFunction.boundaryFieldRef() = 0.0;
401 
402  return tstreamFunction;
403 }
404 
405 
406 bool Foam::functionObjects::streamFunction::calc()
407 {
408  const auto* phiPtr = findObject<surfaceScalarField>(fieldName_);
409 
410  if (phiPtr)
411  {
412  const surfaceScalarField& phi = *phiPtr;
413 
414  return store(resultName_, calc(phi));
415  }
416 
417  return false;
418 }
419 
420 
421 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
422 
424 (
425  const word& name,
426  const Time& runTime,
427  const dictionary& dict
428 )
429 :
430  fieldExpression(name, runTime, dict, "phi")
431 {
432  setResultName(typeName, "phi");
433 
434  const label nD = mesh_.nGeometricD();
435 
436  if (nD != 2)
437  {
439  << "Case is not 2D, stream-function cannot be computed"
440  << exit(FatalError);
441  }
442 }
443 
444 
445 // ************************************************************************* //
Foam::surfaceFields.
dimensionedScalar sign(const dimensionedScalar &ds)
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
streamFunction(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
uint8_t direction
Definition: direction.H:46
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const cellList & cells() const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:876
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
UList< face > faceUList
UList of faces.
Definition: faceListFwd.H:43
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const word & type() const =0
Runtime type information.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
vector point
Point is a vector.
Definition: point.H:37
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
const vectorField & faceAreas() const
const dimensionedScalar c
Speed of light in a vacuum.
const polyBoundaryMesh & patches
#define Log
Definition: PDRblock.C:28
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:192
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool found
const fvMesh & mesh_
Reference to the fvMesh.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const Time & time_
Reference to the time database.