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 | Copyright (C) 2011-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 Application
25  streamFunction
26 
27 Description
28  Calculates and writes the stream function of velocity field U at each
29  time.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "pointFields.H"
35 #include "emptyPolyPatch.H"
36 #include "symmetryPlanePolyPatch.H"
37 #include "symmetryPolyPatch.H"
38 #include "wedgePolyPatch.H"
39 #include "OSspecific.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  timeSelector::addOptions();
46  #include "addRegionOption.H"
47 
48  #include "setRootCase.H"
49  #include "createTime.H"
50 
51  instantList timeDirs = timeSelector::select0(runTime, args);
52 
53  #include "createNamedMesh.H"
54 
55  label nD = mesh.nGeometricD();
56 
57  if (nD != 2)
58  {
60  << "Case is not 2D, stream-function cannot be computed"
61  << exit(FatalError);
62  }
63 
64  Vector<label> slabNormal((Vector<label>::one - mesh.geometricD())/2);
65  const direction slabDir
66  (
67  slabNormal
69  );
70 
71  scalar thickness = vector(slabNormal) & mesh.bounds().span();
72 
73  const pointMesh& pMesh = pointMesh::New(mesh);
74 
75  forAll(timeDirs, timeI)
76  {
77  runTime.setTime(timeDirs[timeI], timeI);
78 
79  Info<< nl << "Time: " << runTime.timeName() << endl;
80 
81  IOobject phiHeader
82  (
83  "phi",
84  runTime.timeName(),
85  mesh,
86  IOobject::NO_READ
87  );
88 
89  if (phiHeader.headerOk())
90  {
91  mesh.readUpdate();
92 
93  Info<< nl << "Reading field phi" << endl;
94 
96  (
97  IOobject
98  (
99  "phi",
100  runTime.timeName(),
101  mesh,
102  IOobject::MUST_READ,
103  IOobject::NO_WRITE
104  ),
105  mesh
106  );
107 
108  pointScalarField streamFunction
109  (
110  IOobject
111  (
112  "streamFunction",
113  runTime.timeName(),
114  mesh,
115  IOobject::NO_READ,
116  IOobject::NO_WRITE
117  ),
118  pMesh,
119  dimensionedScalar("zero", phi.dimensions(), 0.0)
120  );
121 
122  labelList visitedPoint(mesh.nPoints());
123  forAll(visitedPoint, pointi)
124  {
125  visitedPoint[pointi] = 0;
126  }
127  label nVisited = 0;
128  label nVisitedOld = 0;
129 
130  const faceUList& faces = mesh.faces();
131  const pointField& points = mesh.points();
132 
133  label nInternalFaces = mesh.nInternalFaces();
134 
135  vectorField unitAreas(mesh.faceAreas());
136  unitAreas /= mag(unitAreas);
137 
138  const polyPatchList& patches = mesh.boundaryMesh();
139 
140  bool finished = true;
141 
142  // Find the boundary face with zero flux. set the stream function
143  // to zero on that face
144  bool found = false;
145 
146  do
147  {
148  found = false;
149 
150  forAll(patches, patchi)
151  {
152  const primitivePatch& bouFaces = patches[patchi];
153 
154  if (!isType<emptyPolyPatch>(patches[patchi]))
155  {
156  forAll(bouFaces, facei)
157  {
158  if
159  (
160  magSqr(phi.boundaryField()[patchi][facei])
161  < SMALL
162  )
163  {
164  const labelList& zeroPoints = bouFaces[facei];
165 
166  // Zero flux face found
167  found = true;
168 
169  forAll(zeroPoints, pointi)
170  {
171  if (visitedPoint[zeroPoints[pointi]] == 1)
172  {
173  found = false;
174  break;
175  }
176  }
177 
178  if (found)
179  {
180  Info<< "Zero face: patch: " << patchi
181  << " face: " << facei << endl;
182 
183  forAll(zeroPoints, pointi)
184  {
185  streamFunction[zeroPoints[pointi]] = 0;
186  visitedPoint[zeroPoints[pointi]] = 1;
187  nVisited++;
188  }
189 
190  break;
191  }
192  }
193  }
194  }
195 
196  if (found) break;
197  }
198 
199  if (!found)
200  {
201  Info<< "zero flux boundary face not found. "
202  << "Using cell as a reference."
203  << endl;
204 
205  const cellList& c = mesh.cells();
206 
207  forAll(c, cI)
208  {
209  labelList zeroPoints = c[cI].labels(mesh.faces());
210 
211  bool found = true;
212 
213  forAll(zeroPoints, pointi)
214  {
215  if (visitedPoint[zeroPoints[pointi]] == 1)
216  {
217  found = false;
218  break;
219  }
220  }
221 
222  if (found)
223  {
224  forAll(zeroPoints, pointi)
225  {
226  streamFunction[zeroPoints[pointi]] = 0.0;
227  visitedPoint[zeroPoints[pointi]] = 1;
228  nVisited++;
229  }
230 
231  break;
232  }
233  else
234  {
236  << "Cannot find initialisation face or a cell."
237  << abort(FatalError);
238  }
239  }
240  }
241 
242  // Loop through all faces. If one of the points on
243  // the face has the streamfunction value different
244  // from -1, all points with -1 ont that face have the
245  // streamfunction value equal to the face flux in
246  // that point plus the value in the visited point
247  do
248  {
249  finished = true;
250 
251  for
252  (
253  label facei = nInternalFaces;
254  facei<faces.size();
255  facei++
256  )
257  {
258  const labelList& curBPoints = faces[facei];
259  bool bPointFound = false;
260 
261  scalar currentBStream = 0.0;
262  vector currentBStreamPoint(0, 0, 0);
263 
264  forAll(curBPoints, pointi)
265  {
266  // Check if the point has been visited
267  if (visitedPoint[curBPoints[pointi]] == 1)
268  {
269  // The point has been visited
270  currentBStream =
271  streamFunction[curBPoints[pointi]];
272  currentBStreamPoint =
273  points[curBPoints[pointi]];
274 
275  bPointFound = true;
276 
277  break;
278  }
279  }
280 
281  if (bPointFound)
282  {
283  // Sort out other points on the face
284  forAll(curBPoints, pointi)
285  {
286  // Check if the point has been visited
287  if (visitedPoint[curBPoints[pointi]] == 0)
288  {
289  label patchNo =
290  mesh.boundaryMesh().whichPatch(facei);
291 
292  if
293  (
294  !isType<emptyPolyPatch>
295  (patches[patchNo])
296  && !isType<symmetryPlanePolyPatch>
297  (patches[patchNo])
298  && !isType<symmetryPolyPatch>
299  (patches[patchNo])
300  && !isType<wedgePolyPatch>
301  (patches[patchNo])
302  )
303  {
304  label faceNo =
305  mesh.boundaryMesh()[patchNo]
306  .whichFace(facei);
307 
308  vector edgeHat =
309  points[curBPoints[pointi]]
310  - currentBStreamPoint;
311  edgeHat.replace(slabDir, 0);
312  edgeHat /= mag(edgeHat);
313 
314  vector nHat = unitAreas[facei];
315 
316  if (edgeHat.y() > VSMALL)
317  {
318  visitedPoint[curBPoints[pointi]] =
319  1;
320  nVisited++;
321 
322  streamFunction[curBPoints[pointi]]
323  =
324  currentBStream
325  + phi.boundaryField()
326  [patchNo][faceNo]
327  *sign(nHat.x());
328  }
329  else if (edgeHat.y() < -VSMALL)
330  {
331  visitedPoint[curBPoints[pointi]] =
332  1;
333  nVisited++;
334 
335  streamFunction[curBPoints[pointi]]
336  =
337  currentBStream
338  - phi.boundaryField()
339  [patchNo][faceNo]
340  *sign(nHat.x());
341  }
342  else
343  {
344  if (edgeHat.x() > VSMALL)
345  {
346  visitedPoint
347  [curBPoints[pointi]] = 1;
348  nVisited++;
349 
350  streamFunction
351  [curBPoints[pointi]] =
352  currentBStream
353  + phi.boundaryField()
354  [patchNo][faceNo]
355  *sign(nHat.y());
356  }
357  else if (edgeHat.x() < -VSMALL)
358  {
359  visitedPoint
360  [curBPoints[pointi]] = 1;
361  nVisited++;
362 
363  streamFunction
364  [curBPoints[pointi]] =
365  currentBStream
366  - phi.boundaryField()
367  [patchNo][faceNo]
368  *sign(nHat.y());
369  }
370  }
371  }
372  }
373  }
374  }
375  else
376  {
377  finished = false;
378  }
379  }
380 
381  for (label facei=0; facei<nInternalFaces; facei++)
382  {
383  // Get the list of point labels for the face
384  const labelList& curPoints = faces[facei];
385 
386  bool pointFound = false;
387 
388  scalar currentStream = 0.0;
389  point currentStreamPoint(0, 0, 0);
390 
391  forAll(curPoints, pointi)
392  {
393  // Check if the point has been visited
394  if (visitedPoint[curPoints[pointi]] == 1)
395  {
396  // The point has been visited
397  currentStream =
398  streamFunction[curPoints[pointi]];
399  currentStreamPoint =
400  points[curPoints[pointi]];
401  pointFound = true;
402 
403  break;
404  }
405  }
406 
407  if (pointFound)
408  {
409  // Sort out other points on the face
410  forAll(curPoints, pointi)
411  {
412  // Check if the point has been visited
413  if (visitedPoint[curPoints[pointi]] == 0)
414  {
415  vector edgeHat =
416  points[curPoints[pointi]]
417  - currentStreamPoint;
418 
419  edgeHat.replace(slabDir, 0);
420  edgeHat /= mag(edgeHat);
421 
422  vector nHat = unitAreas[facei];
423 
424  if (edgeHat.y() > VSMALL)
425  {
426  visitedPoint[curPoints[pointi]] = 1;
427  nVisited++;
428 
429  streamFunction[curPoints[pointi]] =
430  currentStream
431  + phi[facei]*sign(nHat.x());
432  }
433  else if (edgeHat.y() < -VSMALL)
434  {
435  visitedPoint[curPoints[pointi]] = 1;
436  nVisited++;
437 
438  streamFunction[curPoints[pointi]] =
439  currentStream
440  - phi[facei]*sign(nHat.x());
441  }
442  }
443  }
444  }
445  else
446  {
447  finished = false;
448  }
449  }
450 
451  Info<< ".";
452 
453  if (nVisited == nVisitedOld)
454  {
455  // Find new seed. This must be a
456  // multiply connected domain
457  Info<< nl << "Exhausted a seed. Looking for new seed "
458  << "(this is correct for multiply connected "
459  << "domains).";
460 
461  break;
462  }
463  else
464  {
465  nVisitedOld = nVisited;
466  }
467  } while (!finished);
468 
469  Info<< endl;
470  } while (!finished);
471 
472  // Normalise the stream-function by the 2D mesh thickness
473  streamFunction /= thickness;
474  streamFunction.boundaryFieldRef() = 0.0;
475  streamFunction.write();
476  }
477  else
478  {
480  << "Flux field does not exist."
481  << " Stream function not calculated" << endl;
482  }
483  }
484 
485  Info<< "\nEnd\n" << endl;
486 
487  return 0;
488 }
489 
490 
491 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
dimensionedScalar sign(const dimensionedScalar &ds)
surfaceScalarField & phi
List< instant > instantList
List of instants.
Definition: instantList.H:42
#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
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
volVectorField vectorField(fieldObject, mesh)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void replace(const direction, const Cmpt &)
Definition: VectorSpaceI.H:159
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
UList< face > faceUList
Definition: faceListFwd.H:41
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:262
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
label patchi
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Foam::argList args(argc, argv)
List< cell > cellList
list of cells
Definition: cellList.H:42