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